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 MatHYPRE_IJMatrixFastCopy(Mat,HYPRE_IJMatrix); 20 EXTERN PetscErrorCode VecHYPRE_IJVectorCreate(Vec,HYPRE_IJVector*); 21 22 /* 23 Private context (data structure) for the preconditioner. 24 */ 25 typedef struct { 26 HYPRE_Solver hsolver; 27 HYPRE_IJMatrix ij; 28 HYPRE_IJVector b,x; 29 30 PetscErrorCode (*destroy)(HYPRE_Solver); 31 PetscErrorCode (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector); 32 PetscErrorCode (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector); 33 34 MPI_Comm comm_hypre; 35 char *hypre_type; 36 37 /* options for Pilut and BoomerAMG*/ 38 int maxiter; 39 double tol; 40 41 /* options for Pilut */ 42 int factorrowsize; 43 44 /* options for ParaSails */ 45 int nlevels; 46 double threshhold; 47 double filter; 48 int sym; 49 double loadbal; 50 int logging; 51 int ruse; 52 int symt; 53 54 /* options for Euclid */ 55 PetscTruth bjilu; 56 int levels; 57 58 /* options for Euclid and BoomerAMG */ 59 PetscTruth printstatistics; 60 61 /* options for BoomerAMG */ 62 int cycletype; 63 int maxlevels; 64 double strongthreshold; 65 double maxrowsum; 66 int gridsweeps[3]; 67 int coarsentype; 68 int measuretype; 69 int relaxtype[3]; 70 double relaxweight; 71 double outerrelaxweight; 72 int relaxorder; 73 double truncfactor; 74 PetscTruth applyrichardson; 75 int pmax; 76 int interptype; 77 int agg_nl; 78 int agg_num_paths; 79 int nodal_coarsen; 80 PetscTruth nodal_relax; 81 int nodal_relax_levels; 82 } PC_HYPRE; 83 84 85 #undef __FUNCT__ 86 #define __FUNCT__ "PCSetUp_HYPRE" 87 static PetscErrorCode PCSetUp_HYPRE(PC pc) 88 { 89 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 90 PetscErrorCode ierr; 91 HYPRE_ParCSRMatrix hmat; 92 HYPRE_ParVector bv,xv; 93 PetscInt bs; 94 int hierr; 95 96 PetscFunctionBegin; 97 if (!jac->hypre_type) { 98 ierr = PCHYPRESetType(pc,"boomeramg");CHKERRQ(ierr); 99 } 100 101 if (pc->setupcalled) { 102 /* always destroy the old matrix and create a new memory; 103 hope this does not churn the memory too much. The problem 104 is I do not know if it is possible to put the matrix back to 105 its initial state so that we can directly copy the values 106 the second time through. */ 107 ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); 108 jac->ij = 0; 109 } 110 111 if (!jac->ij) { /* create the matrix the first time through */ 112 ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr); 113 } 114 if (!jac->b) { /* create the vectors the first time through */ 115 Vec x,b; 116 ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr); 117 ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr); 118 ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr); 119 ierr = VecDestroy(x);CHKERRQ(ierr); 120 ierr = VecDestroy(b);CHKERRQ(ierr); 121 } 122 123 /* special case for BoomerAMG */ 124 if (jac->setup == HYPRE_BoomerAMGSetup) { 125 ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr); 126 if (bs > 1) { 127 ierr = HYPRE_BoomerAMGSetNumFunctions(jac->hsolver,bs);CHKERRQ(ierr); 128 } 129 }; 130 ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr); 131 ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 132 ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&bv);CHKERRQ(ierr); 133 ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&xv);CHKERRQ(ierr); 134 hierr = (*jac->setup)(jac->hsolver,hmat,bv,xv); 135 if (hierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE setup, error code %d",hierr); 136 PetscFunctionReturn(0); 137 } 138 139 /* 140 Replaces the address where the HYPRE vector points to its data with the address of 141 PETSc's data. Saves the old address so it can be reset when we are finished with it. 142 Allows use to get the data into a HYPRE vector without the cost of memcopies 143 */ 144 #define HYPREReplacePointer(b,newvalue,savedvalue) {\ 145 hypre_ParVector *par_vector = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\ 146 hypre_Vector *local_vector = hypre_ParVectorLocalVector(par_vector);\ 147 savedvalue = local_vector->data;\ 148 local_vector->data = newvalue;} 149 150 #undef __FUNCT__ 151 #define __FUNCT__ "PCApply_HYPRE" 152 static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x) 153 { 154 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 155 PetscErrorCode ierr; 156 HYPRE_ParCSRMatrix hmat; 157 PetscScalar *bv,*xv; 158 HYPRE_ParVector jbv,jxv; 159 PetscScalar *sbv,*sxv; 160 int hierr; 161 162 PetscFunctionBegin; 163 if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);} 164 ierr = VecGetArray(b,&bv);CHKERRQ(ierr); 165 ierr = VecGetArray(x,&xv);CHKERRQ(ierr); 166 HYPREReplacePointer(jac->b,bv,sbv); 167 HYPREReplacePointer(jac->x,xv,sxv); 168 169 ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 170 ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr); 171 ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr); 172 hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv); 173 174 /*if (hierr && (hierr != HYPRE_ERROR_CONV || jac->solve != HYPRE_BoomerAMGSolve))SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 175 */ 176 /* error code of HYPRE_ERROR_CONV means convergence not achieved - if 177 the tolerance is set to 0.0 (the default), a convergence error will 178 not occur (so we may not want to overide the conv. error here?*/ 179 if (hierr && hierr != HYPRE_ERROR_CONV) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 180 if (hierr) hypre__global_error = 0; 181 182 183 HYPREReplacePointer(jac->b,sbv,bv); 184 HYPREReplacePointer(jac->x,sxv,xv); 185 ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr); 186 ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr); 187 PetscFunctionReturn(0); 188 } 189 190 #undef __FUNCT__ 191 #define __FUNCT__ "PCDestroy_HYPRE" 192 static PetscErrorCode PCDestroy_HYPRE(PC pc) 193 { 194 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 195 PetscErrorCode ierr; 196 197 PetscFunctionBegin; 198 if (jac->ij) { ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); } 199 if (jac->b) { ierr = HYPRE_IJVectorDestroy(jac->b);CHKERRQ(ierr); } 200 if (jac->x) { ierr = HYPRE_IJVectorDestroy(jac->x);CHKERRQ(ierr); } 201 if (jac->destroy) { ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr); } 202 ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr); 203 if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);} 204 ierr = PetscFree(jac);CHKERRQ(ierr); 205 206 ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr); 207 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr); 208 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr); 209 PetscFunctionReturn(0); 210 } 211 212 /* --------------------------------------------------------------------------------------------*/ 213 #undef __FUNCT__ 214 #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut" 215 static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc) 216 { 217 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 218 PetscErrorCode ierr; 219 PetscTruth flag; 220 221 PetscFunctionBegin; 222 ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr); 223 ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr); 224 if (flag) { 225 ierr = HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 226 } 227 ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr); 228 if (flag) { 229 ierr = HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol);CHKERRQ(ierr); 230 } 231 ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr); 232 if (flag) { 233 ierr = HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize);CHKERRQ(ierr); 234 } 235 ierr = PetscOptionsTail();CHKERRQ(ierr); 236 PetscFunctionReturn(0); 237 } 238 239 #undef __FUNCT__ 240 #define __FUNCT__ "PCView_HYPRE_Pilut" 241 static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer) 242 { 243 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 244 PetscErrorCode ierr; 245 PetscTruth iascii; 246 247 PetscFunctionBegin; 248 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 249 if (iascii) { 250 ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut preconditioning\n");CHKERRQ(ierr); 251 if (jac->maxiter != PETSC_DEFAULT) { 252 ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr); 253 } else { 254 ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr); 255 } 256 if (jac->tol != PETSC_DEFAULT) { 257 ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr); 258 } else { 259 ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr); 260 } 261 if (jac->factorrowsize != PETSC_DEFAULT) { 262 ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr); 263 } else { 264 ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default factor row size \n");CHKERRQ(ierr); 265 } 266 } 267 PetscFunctionReturn(0); 268 } 269 270 /* --------------------------------------------------------------------------------------------*/ 271 #undef __FUNCT__ 272 #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid" 273 static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc) 274 { 275 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 276 PetscErrorCode ierr; 277 PetscTruth flag; 278 char *args[8],levels[16]; 279 PetscInt cnt = 0; 280 281 PetscFunctionBegin; 282 ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr); 283 ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr); 284 if (flag) { 285 if (jac->levels < 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels); 286 sprintf(levels,"%d",jac->levels); 287 args[cnt++] = (char*)"-level"; args[cnt++] = levels; 288 } 289 ierr = PetscOptionsTruth("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr); 290 if (jac->bjilu) { 291 args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1"; 292 } 293 294 ierr = PetscOptionsTruth("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr); 295 if (jac->printstatistics) { 296 args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1"; 297 args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1"; 298 } 299 ierr = PetscOptionsTail();CHKERRQ(ierr); 300 if (cnt) { 301 ierr = HYPRE_EuclidSetParams(jac->hsolver,cnt,args);CHKERRQ(ierr); 302 } 303 PetscFunctionReturn(0); 304 } 305 306 #undef __FUNCT__ 307 #define __FUNCT__ "PCView_HYPRE_Euclid" 308 static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer) 309 { 310 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 311 PetscErrorCode ierr; 312 PetscTruth iascii; 313 314 PetscFunctionBegin; 315 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 316 if (iascii) { 317 ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid preconditioning\n");CHKERRQ(ierr); 318 ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr); 319 if (jac->bjilu) { 320 ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr); 321 } 322 } 323 PetscFunctionReturn(0); 324 } 325 326 /* --------------------------------------------------------------------------------------------*/ 327 328 #undef __FUNCT__ 329 #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG" 330 static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x) 331 { 332 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 333 PetscErrorCode ierr; 334 HYPRE_ParCSRMatrix hmat; 335 PetscScalar *bv,*xv; 336 HYPRE_ParVector jbv,jxv; 337 PetscScalar *sbv,*sxv; 338 int hierr; 339 340 PetscFunctionBegin; 341 ierr = VecSet(x,0.0);CHKERRQ(ierr); 342 ierr = VecGetArray(b,&bv);CHKERRQ(ierr); 343 ierr = VecGetArray(x,&xv);CHKERRQ(ierr); 344 HYPREReplacePointer(jac->b,bv,sbv); 345 HYPREReplacePointer(jac->x,xv,sxv); 346 347 ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 348 ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr); 349 ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr); 350 351 hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv); 352 /* error code of 1 in BoomerAMG merely means convergence not achieved */ 353 if (hierr && (hierr != 1)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 354 if (hierr) hypre__global_error = 0; 355 356 HYPREReplacePointer(jac->b,sbv,bv); 357 HYPREReplacePointer(jac->x,sxv,xv); 358 ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr); 359 ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr); 360 PetscFunctionReturn(0); 361 } 362 363 static const char *HYPREBoomerAMGCycleType[] = {"","V","W"}; 364 static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"}; 365 static const char *HYPREBoomerAMGMeasureType[] = {"local","global"}; 366 static const char *HYPREBoomerAMGRelaxType[] = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi", 367 "","","Gaussian-elimination"}; 368 static const char *HYPREBoomerAMGInterpType[] = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i", 369 "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"}; 370 #undef __FUNCT__ 371 #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG" 372 static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc) 373 { 374 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 375 PetscErrorCode ierr; 376 int n,indx; 377 PetscTruth flg, tmp_truth; 378 double tmpdbl, twodbl[2]; 379 380 PetscFunctionBegin; 381 ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr); 382 ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr); 383 if (flg) { 384 jac->cycletype = indx; 385 ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr); 386 } 387 ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr); 388 if (flg) { 389 if (jac->maxlevels < 2) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels); 390 ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr); 391 } 392 ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr); 393 if (flg) { 394 if (jac->maxiter < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter); 395 ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 396 } 397 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); 398 if (flg) { 399 if (jac->tol < 0.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be greater than or equal to zero",jac->tol); 400 ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 401 } 402 403 ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr); 404 if (flg) { 405 if (jac->truncfactor < 0.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor); 406 ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr); 407 } 408 409 ierr = PetscOptionsInt("-pc_hypre_boomeramg_P_max","Max elements per row for interpolation operator ( 0=unlimited )","None",jac->pmax,&jac->pmax,&flg);CHKERRQ(ierr); 410 if (flg) { 411 if (jac->pmax < 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"P_max %G must be greater than or equal to zero",jac->pmax); 412 ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr); 413 } 414 415 ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr); 416 if (flg) { 417 if (jac->agg_nl < 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %G must be greater than or equal to zero",jac->agg_nl); 418 419 ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl);CHKERRQ(ierr); 420 } 421 422 423 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); 424 if (flg) { 425 if (jac->agg_num_paths < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of paths %G must be greater than or equal to 1",jac->agg_num_paths); 426 427 ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr); 428 } 429 430 431 ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr); 432 if (flg) { 433 if (jac->strongthreshold < 0.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold); 434 ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr); 435 } 436 ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr); 437 if (flg) { 438 if (jac->maxrowsum < 0.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum); 439 if (jac->maxrowsum > 1.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum); 440 ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr); 441 } 442 443 /* Grid sweeps */ 444 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); 445 if (flg) { 446 ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver,indx);CHKERRQ(ierr); 447 /* modify the jac structure so we can view the updated options with PC_View */ 448 jac->gridsweeps[0] = indx; 449 jac->gridsweeps[1] = indx; 450 /*defaults coarse to 1 */ 451 jac->gridsweeps[2] = 1; 452 } 453 454 ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx ,&flg);CHKERRQ(ierr); 455 if (flg) { 456 ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 1);CHKERRQ(ierr); 457 jac->gridsweeps[0] = indx; 458 } 459 ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr); 460 if (flg) { 461 ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 2);CHKERRQ(ierr); 462 jac->gridsweeps[1] = indx; 463 } 464 ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr); 465 if (flg) { 466 ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 3);CHKERRQ(ierr); 467 jac->gridsweeps[2] = indx; 468 } 469 470 /* Relax type */ 471 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); 472 if (flg) { 473 jac->relaxtype[0] = jac->relaxtype[1] = indx; 474 ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, indx);CHKERRQ(ierr); 475 /* by default, coarse type set to 9 */ 476 jac->relaxtype[2] = 9; 477 478 } 479 ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr); 480 if (flg) { 481 jac->relaxtype[0] = indx; 482 ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 1);CHKERRQ(ierr); 483 } 484 ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr); 485 if (flg) { 486 jac->relaxtype[1] = indx; 487 ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 2);CHKERRQ(ierr); 488 } 489 ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr); 490 if (flg) { 491 jac->relaxtype[2] = indx; 492 ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 3);CHKERRQ(ierr); 493 } 494 495 /* Relaxation Weight */ 496 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); 497 if (flg) { 498 ierr = HYPRE_BoomerAMGSetRelaxWt(jac->hsolver,tmpdbl);CHKERRQ(ierr); 499 jac->relaxweight = tmpdbl; 500 } 501 502 n=2; 503 twodbl[0] = twodbl[1] = 1.0; 504 ierr = PetscOptionsRealArray("-pc_hypre_boomeramg_relax_weight_level","Set the relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);CHKERRQ(ierr); 505 if (flg) { 506 if (n == 2) { 507 indx = (int)PetscAbsReal(twodbl[1]); 508 ierr = HYPRE_BoomerAMGSetLevelRelaxWt(jac->hsolver,twodbl[0],indx);CHKERRQ(ierr); 509 } else { 510 SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relax weight level: you must provide 2 values separated by a comma (and no space), you provided %d",n); 511 } 512 } 513 514 /* Outer relaxation Weight */ 515 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); 516 if (flg) { 517 ierr = HYPRE_BoomerAMGSetOuterWt( jac->hsolver, tmpdbl);CHKERRQ(ierr); 518 jac->outerrelaxweight = tmpdbl; 519 } 520 521 n=2; 522 twodbl[0] = twodbl[1] = 1.0; 523 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); 524 if (flg) { 525 if (n == 2) { 526 indx = (int)PetscAbsReal(twodbl[1]); 527 ierr = HYPRE_BoomerAMGSetLevelOuterWt( jac->hsolver, twodbl[0], indx);CHKERRQ(ierr); 528 } else { 529 SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relax weight outer level: You must provide 2 values separated by a comma (and no space), you provided %d",n); 530 } 531 } 532 533 /* the Relax Order */ 534 ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 535 536 if (flg) { 537 jac->relaxorder = 0; 538 ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr); 539 } 540 ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr); 541 if (flg) { 542 jac->measuretype = indx; 543 ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr); 544 } 545 /* update list length 3/07 */ 546 ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,11,HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr); 547 if (flg) { 548 jac->coarsentype = indx; 549 ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr); 550 } 551 552 /* new 3/07 */ 553 ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,14,HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr); 554 if (flg) { 555 jac->interptype = indx; 556 ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr); 557 } 558 559 flg = PETSC_FALSE; 560 ierr = PetscOptionsTruth("-pc_hypre_boomeramg_print_statistics","Print statistics","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 561 if (flg) { 562 int level=3; 563 jac->printstatistics = PETSC_TRUE; 564 ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr); 565 ierr = HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,level);CHKERRQ(ierr); 566 } 567 568 flg = PETSC_FALSE; 569 ierr = PetscOptionsTruth("-pc_hypre_boomeramg_print_debug","Print debug information","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 570 if (flg) { 571 int level=3; 572 jac->printstatistics = PETSC_TRUE; 573 ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,PETSC_NULL);CHKERRQ(ierr); 574 ierr = HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,level);CHKERRQ(ierr); 575 } 576 577 ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 578 if (flg && tmp_truth) { 579 jac->nodal_coarsen = 1; 580 ierr = HYPRE_BoomerAMGSetNodal(jac->hsolver,1);CHKERRQ(ierr); 581 } 582 583 ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 584 if (flg && tmp_truth) { 585 PetscInt tmp_int; 586 ierr = PetscOptionsInt( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);CHKERRQ(ierr); 587 if (flg) jac->nodal_relax_levels = tmp_int; 588 ierr = HYPRE_BoomerAMGSetSmoothType(jac->hsolver,6);CHKERRQ(ierr); 589 ierr = HYPRE_BoomerAMGSetDomainType(jac->hsolver,1);CHKERRQ(ierr); 590 ierr = HYPRE_BoomerAMGSetOverlap(jac->hsolver,0);CHKERRQ(ierr); 591 ierr = HYPRE_BoomerAMGSetSmoothNumLevels(jac->hsolver,jac->nodal_relax_levels);CHKERRQ(ierr); 592 } 593 594 ierr = PetscOptionsTail();CHKERRQ(ierr); 595 PetscFunctionReturn(0); 596 } 597 598 #undef __FUNCT__ 599 #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG" 600 static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason) 601 { 602 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 603 PetscErrorCode ierr; 604 int oits; 605 606 PetscFunctionBegin; 607 ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its*jac->maxiter);CHKERRQ(ierr); 608 ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,rtol);CHKERRQ(ierr); 609 jac->applyrichardson = PETSC_TRUE; 610 ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr); 611 jac->applyrichardson = PETSC_FALSE; 612 ierr = HYPRE_BoomerAMGGetNumIterations(jac->hsolver,&oits);CHKERRQ(ierr); 613 *outits = oits; 614 if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS; 615 else *reason = PCRICHARDSON_CONVERGED_RTOL; 616 ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 617 ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 618 PetscFunctionReturn(0); 619 } 620 621 622 #undef __FUNCT__ 623 #define __FUNCT__ "PCView_HYPRE_BoomerAMG" 624 static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer) 625 { 626 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 627 PetscErrorCode ierr; 628 PetscTruth iascii; 629 630 PetscFunctionBegin; 631 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 632 if (iascii) { 633 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr); 634 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr); 635 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr); 636 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr); 637 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr); 638 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr); 639 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);CHKERRQ(ierr); 640 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr); 641 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr); 642 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr); 643 644 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr); 645 646 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps down %d\n",jac->gridsweeps[0]);CHKERRQ(ierr); 647 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps up %d\n",jac->gridsweeps[1]);CHKERRQ(ierr); 648 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps on coarse %d\n",jac->gridsweeps[2]);CHKERRQ(ierr); 649 650 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax down %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr); 651 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax up %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr); 652 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax on coarse %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr); 653 654 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax weight (all) %G\n",jac->relaxweight);CHKERRQ(ierr); 655 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr); 656 657 if (jac->relaxorder) { 658 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr); 659 } else { 660 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr); 661 } 662 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Measure type %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr); 663 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Coarsen type %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr); 664 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation type %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr); 665 if (jac->nodal_coarsen) { 666 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");CHKERRQ(ierr); 667 } 668 if (jac->nodal_relax) { 669 ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);CHKERRQ(ierr); 670 } 671 } 672 PetscFunctionReturn(0); 673 } 674 675 /* --------------------------------------------------------------------------------------------*/ 676 #undef __FUNCT__ 677 #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails" 678 static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc) 679 { 680 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 681 PetscErrorCode ierr; 682 int indx; 683 PetscTruth flag; 684 const char *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"}; 685 686 PetscFunctionBegin; 687 ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr); 688 ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr); 689 ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr); 690 if (flag) { 691 ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr); 692 } 693 694 ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr); 695 if (flag) { 696 ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr); 697 } 698 699 ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr); 700 if (flag) { 701 ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr); 702 } 703 704 ierr = PetscOptionsTruth("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,&flag);CHKERRQ(ierr); 705 if (flag) { 706 ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr); 707 } 708 709 ierr = PetscOptionsTruth("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,&flag);CHKERRQ(ierr); 710 if (flag) { 711 ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr); 712 } 713 714 ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);CHKERRQ(ierr); 715 if (flag) { 716 jac->symt = indx; 717 ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr); 718 } 719 720 ierr = PetscOptionsTail();CHKERRQ(ierr); 721 PetscFunctionReturn(0); 722 } 723 724 #undef __FUNCT__ 725 #define __FUNCT__ "PCView_HYPRE_ParaSails" 726 static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer) 727 { 728 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 729 PetscErrorCode ierr; 730 PetscTruth iascii; 731 const char *symt = 0;; 732 733 PetscFunctionBegin; 734 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 735 if (iascii) { 736 ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails preconditioning\n");CHKERRQ(ierr); 737 ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr); 738 ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr); 739 ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr); 740 ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr); 741 ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: reuse nonzero structure %s\n",PetscTruths[jac->ruse]);CHKERRQ(ierr); 742 ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: print info to screen %s\n",PetscTruths[jac->logging]);CHKERRQ(ierr); 743 if (!jac->symt) { 744 symt = "nonsymmetric matrix and preconditioner"; 745 } else if (jac->symt == 1) { 746 symt = "SPD matrix and preconditioner"; 747 } else if (jac->symt == 2) { 748 symt = "nonsymmetric matrix but SPD preconditioner"; 749 } else { 750 SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt); 751 } 752 ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr); 753 } 754 PetscFunctionReturn(0); 755 } 756 /* ---------------------------------------------------------------------------------*/ 757 758 EXTERN_C_BEGIN 759 #undef __FUNCT__ 760 #define __FUNCT__ "PCHYPREGetType_HYPRE" 761 PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType_HYPRE(PC pc,const char *name[]) 762 { 763 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 764 765 PetscFunctionBegin; 766 *name = jac->hypre_type; 767 PetscFunctionReturn(0); 768 } 769 EXTERN_C_END 770 771 EXTERN_C_BEGIN 772 #undef __FUNCT__ 773 #define __FUNCT__ "PCHYPRESetType_HYPRE" 774 PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType_HYPRE(PC pc,const char name[]) 775 { 776 PC_HYPRE *jac = (PC_HYPRE*)pc->data; 777 PetscErrorCode ierr; 778 PetscTruth flag; 779 780 PetscFunctionBegin; 781 if (jac->hypre_type) { 782 ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr); 783 if (!flag) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set"); 784 PetscFunctionReturn(0); 785 } else { 786 ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr); 787 } 788 789 jac->maxiter = PETSC_DEFAULT; 790 jac->tol = PETSC_DEFAULT; 791 jac->printstatistics = PetscLogPrintInfo; 792 793 ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr); 794 if (flag) { 795 ierr = HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver); 796 pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut; 797 pc->ops->view = PCView_HYPRE_Pilut; 798 jac->destroy = HYPRE_ParCSRPilutDestroy; 799 jac->setup = HYPRE_ParCSRPilutSetup; 800 jac->solve = HYPRE_ParCSRPilutSolve; 801 jac->factorrowsize = PETSC_DEFAULT; 802 PetscFunctionReturn(0); 803 } 804 ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr); 805 if (flag) { 806 ierr = HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver); 807 pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails; 808 pc->ops->view = PCView_HYPRE_ParaSails; 809 jac->destroy = HYPRE_ParaSailsDestroy; 810 jac->setup = HYPRE_ParaSailsSetup; 811 jac->solve = HYPRE_ParaSailsSolve; 812 /* initialize */ 813 jac->nlevels = 1; 814 jac->threshhold = .1; 815 jac->filter = .1; 816 jac->loadbal = 0; 817 if (PetscLogPrintInfo) { 818 jac->logging = (int) PETSC_TRUE; 819 } else { 820 jac->logging = (int) PETSC_FALSE; 821 } 822 jac->ruse = (int) PETSC_FALSE; 823 jac->symt = 0; 824 ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr); 825 ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr); 826 ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr); 827 ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr); 828 ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr); 829 ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr); 830 PetscFunctionReturn(0); 831 } 832 ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr); 833 if (flag) { 834 ierr = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver); 835 pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid; 836 pc->ops->view = PCView_HYPRE_Euclid; 837 jac->destroy = HYPRE_EuclidDestroy; 838 jac->setup = HYPRE_EuclidSetup; 839 jac->solve = HYPRE_EuclidSolve; 840 /* initialization */ 841 jac->bjilu = PETSC_FALSE; 842 jac->levels = 1; 843 PetscFunctionReturn(0); 844 } 845 ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr); 846 if (flag) { 847 ierr = HYPRE_BoomerAMGCreate(&jac->hsolver); 848 pc->ops->setfromoptions = PCSetFromOptions_HYPRE_BoomerAMG; 849 pc->ops->view = PCView_HYPRE_BoomerAMG; 850 pc->ops->applytranspose = PCApplyTranspose_HYPRE_BoomerAMG; 851 pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG; 852 jac->destroy = HYPRE_BoomerAMGDestroy; 853 jac->setup = HYPRE_BoomerAMGSetup; 854 jac->solve = HYPRE_BoomerAMGSolve; 855 jac->applyrichardson = PETSC_FALSE; 856 /* these defaults match the hypre defaults */ 857 jac->cycletype = 1; 858 jac->maxlevels = 25; 859 jac->maxiter = 1; 860 jac->tol = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */ 861 jac->truncfactor = 0.0; 862 jac->strongthreshold = .25; 863 jac->maxrowsum = .9; 864 jac->coarsentype = 6; 865 jac->measuretype = 0; 866 jac->gridsweeps[0] = jac->gridsweeps[1] = jac->gridsweeps[2] = 1; 867 jac->relaxtype[0] = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */ 868 jac->relaxtype[2] = 9; /*G.E. */ 869 jac->relaxweight = 1.0; 870 jac->outerrelaxweight = 1.0; 871 jac->relaxorder = 1; 872 jac->interptype = 0; 873 jac->agg_nl = 0; 874 jac->pmax = 0; 875 jac->truncfactor = 0.0; 876 jac->agg_num_paths = 1; 877 878 jac->nodal_coarsen = 0; 879 jac->nodal_relax = PETSC_FALSE; 880 jac->nodal_relax_levels = 1; 881 ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr); 882 ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr); 883 ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 884 ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 885 ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr); 886 ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr); 887 ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr); 888 ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr); 889 ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr); 890 ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr); 891 ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr); 892 ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl); 893 ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr); 894 ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr); 895 ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, jac->relaxtype[0]); /*defaults coarse to 9*/ 896 ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver, jac->gridsweeps[0]); /*defaults coarse to 1 */ 897 PetscFunctionReturn(0); 898 } 899 ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr); 900 jac->hypre_type = PETSC_NULL; 901 SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name); 902 PetscFunctionReturn(0); 903 } 904 EXTERN_C_END 905 906 /* 907 It only gets here if the HYPRE type has not been set before the call to 908 ...SetFromOptions() which actually is most of the time 909 */ 910 #undef __FUNCT__ 911 #define __FUNCT__ "PCSetFromOptions_HYPRE" 912 static PetscErrorCode PCSetFromOptions_HYPRE(PC pc) 913 { 914 PetscErrorCode ierr; 915 int indx; 916 const char *type[] = {"pilut","parasails","boomeramg","euclid"}; 917 PetscTruth flg; 918 919 PetscFunctionBegin; 920 ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr); 921 ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"boomeramg",&indx,&flg);CHKERRQ(ierr); 922 if (flg) { 923 ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr); 924 } else { 925 ierr = PCHYPRESetType_HYPRE(pc,"boomeramg");CHKERRQ(ierr); 926 } 927 if (pc->ops->setfromoptions) { 928 ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr); 929 } 930 ierr = PetscOptionsTail();CHKERRQ(ierr); 931 PetscFunctionReturn(0); 932 } 933 934 #undef __FUNCT__ 935 #define __FUNCT__ "PCHYPRESetType" 936 /*@C 937 PCHYPRESetType - Sets which hypre preconditioner you wish to use 938 939 Input Parameters: 940 + pc - the preconditioner context 941 - name - either pilut, parasails, boomeramg, euclid 942 943 Options Database Keys: 944 -pc_hypre_type - One of pilut, parasails, boomeramg, euclid 945 946 Level: intermediate 947 948 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 949 PCHYPRE 950 951 @*/ 952 PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC pc,const char name[]) 953 { 954 PetscErrorCode ierr,(*f)(PC,const char[]); 955 956 PetscFunctionBegin; 957 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 958 PetscValidCharPointer(name,2); 959 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);CHKERRQ(ierr); 960 if (f) { 961 ierr = (*f)(pc,name);CHKERRQ(ierr); 962 } 963 PetscFunctionReturn(0); 964 } 965 966 #undef __FUNCT__ 967 #define __FUNCT__ "PCHYPREGetType" 968 /*@C 969 PCHYPREGetType - Gets which hypre preconditioner you are using 970 971 Input Parameter: 972 . pc - the preconditioner context 973 974 Output Parameter: 975 . name - either pilut, parasails, boomeramg, euclid 976 977 Level: intermediate 978 979 .seealso: PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC, 980 PCHYPRE 981 982 @*/ 983 PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC pc,const char *name[]) 984 { 985 PetscErrorCode ierr,(*f)(PC,const char*[]); 986 987 PetscFunctionBegin; 988 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 989 PetscValidPointer(name,2); 990 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPREGetType_C",(void (**)(void))&f);CHKERRQ(ierr); 991 if (f) { 992 ierr = (*f)(pc,name);CHKERRQ(ierr); 993 } 994 PetscFunctionReturn(0); 995 } 996 997 /*MC 998 PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre 999 1000 Options Database Keys: 1001 + -pc_hypre_type - One of pilut, parasails, boomeramg, euclid 1002 - Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX 1003 preconditioner 1004 1005 Level: intermediate 1006 1007 Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()), 1008 the many hypre options can ONLY be set via the options database (e.g. the command line 1009 or with PetscOptionsSetValue(), there are no functions to set them) 1010 1011 The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations 1012 (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if 1013 -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner 1014 (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of 1015 iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations 1016 and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10 1017 then AT MOST twenty V-cycles of boomeramg will be called. 1018 1019 Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation 1020 (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry. 1021 Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi. 1022 If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly 1023 and use -ksp_max_it to control the number of V-cycles. 1024 (see the PETSc FAQ.html at the PETSc website under the Documentation tab). 1025 1026 2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option 1027 -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L. 1028 1029 See PCPFMG for access to the hypre Struct PFMG solver 1030 1031 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 1032 PCHYPRESetType(), PCPFMG 1033 1034 M*/ 1035 1036 EXTERN_C_BEGIN 1037 #undef __FUNCT__ 1038 #define __FUNCT__ "PCCreate_HYPRE" 1039 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_HYPRE(PC pc) 1040 { 1041 PC_HYPRE *jac; 1042 PetscErrorCode ierr; 1043 1044 PetscFunctionBegin; 1045 ierr = PetscNewLog(pc,PC_HYPRE,&jac);CHKERRQ(ierr); 1046 pc->data = jac; 1047 pc->ops->destroy = PCDestroy_HYPRE; 1048 pc->ops->setfromoptions = PCSetFromOptions_HYPRE; 1049 pc->ops->setup = PCSetUp_HYPRE; 1050 pc->ops->apply = PCApply_HYPRE; 1051 jac->comm_hypre = MPI_COMM_NULL; 1052 jac->hypre_type = PETSC_NULL; 1053 /* duplicate communicator for hypre */ 1054 ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(jac->comm_hypre));CHKERRQ(ierr); 1055 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr); 1056 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr); 1057 PetscFunctionReturn(0); 1058 } 1059 EXTERN_C_END 1060 1061 /* ---------------------------------------------------------------------------------------------------------------------------------*/ 1062 1063 /* working with a HYPRE_StructMatrix and need access to its data */ 1064 #include "../src/dm/da/utils/mhyp.h" 1065 /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */ 1066 #include "private/matimpl.h" 1067 1068 typedef struct { 1069 MPI_Comm hcomm; /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */ 1070 HYPRE_StructSolver hsolver; 1071 1072 /* keep copy of PFMG options used so may view them */ 1073 int its; 1074 double tol; 1075 int relax_type; 1076 int rap_type; 1077 int num_pre_relax,num_post_relax; 1078 int max_levels; 1079 } PC_PFMG; 1080 1081 #undef __FUNCT__ 1082 #define __FUNCT__ "PCDestroy_PFMG" 1083 PetscErrorCode PCDestroy_PFMG(PC pc) 1084 { 1085 PetscErrorCode ierr; 1086 PC_PFMG *ex = (PC_PFMG*) pc->data; 1087 1088 PetscFunctionBegin; 1089 if (ex->hsolver) {ierr = HYPRE_StructPFMGDestroy(ex->hsolver);CHKERRQ(ierr);} 1090 ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr); 1091 ierr = PetscFree(ex);CHKERRQ(ierr); 1092 PetscFunctionReturn(0); 1093 } 1094 1095 static const char *PFMGRelaxType[] = {"Jacobi","Weighted-Jacobi","symmetric-Red/Black-Gauss-Seidel","Red/Black-Gauss-Seidel"}; 1096 static const char *PFMGRAPType[] = {"Galerkin","non-Galerkin"}; 1097 1098 #undef __FUNCT__ 1099 #define __FUNCT__ "PCView_PFMG" 1100 PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer) 1101 { 1102 PetscErrorCode ierr; 1103 PetscTruth iascii; 1104 PC_PFMG *ex = (PC_PFMG*) pc->data; 1105 1106 PetscFunctionBegin; 1107 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1108 if (iascii) { 1109 ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG preconditioning\n");CHKERRQ(ierr); 1110 ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG: max iterations %d\n",ex->its);CHKERRQ(ierr); 1111 ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr); 1112 ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr); 1113 ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG: RAP type %s\n",PFMGRAPType[ex->rap_type]);CHKERRQ(ierr); 1114 ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr); 1115 ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG: max levels %d\n",ex->max_levels);CHKERRQ(ierr); 1116 } 1117 PetscFunctionReturn(0); 1118 } 1119 1120 1121 #undef __FUNCT__ 1122 #define __FUNCT__ "PCSetFromOptions_PFMG" 1123 PetscErrorCode PCSetFromOptions_PFMG(PC pc) 1124 { 1125 PetscErrorCode ierr; 1126 PC_PFMG *ex = (PC_PFMG*) pc->data; 1127 PetscTruth flg = PETSC_FALSE; 1128 1129 PetscFunctionBegin; 1130 ierr = PetscOptionsHead("PFMG options");CHKERRQ(ierr); 1131 ierr = PetscOptionsTruth("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 1132 if (flg) { 1133 int level=3; 1134 ierr = HYPRE_StructPFMGSetPrintLevel(ex->hsolver,level);CHKERRQ(ierr); 1135 } 1136 ierr = PetscOptionsInt("-pc_pfmg_its","Number of iterations of PFMG to use as preconditioner","HYPRE_StructPFMGSetMaxIter",ex->its,&ex->its,PETSC_NULL);CHKERRQ(ierr); 1137 ierr = HYPRE_StructPFMGSetMaxIter(ex->hsolver,ex->its);CHKERRQ(ierr); 1138 ierr = PetscOptionsInt("-pc_pfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_StructPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,PETSC_NULL);CHKERRQ(ierr); 1139 ierr = HYPRE_StructPFMGSetNumPreRelax(ex->hsolver,ex->num_pre_relax);CHKERRQ(ierr); 1140 ierr = PetscOptionsInt("-pc_pfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_StructPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,PETSC_NULL);CHKERRQ(ierr); 1141 ierr = HYPRE_StructPFMGSetNumPostRelax(ex->hsolver,ex->num_post_relax);CHKERRQ(ierr); 1142 1143 ierr = PetscOptionsInt("-pc_pfmg_max_levels","Max Levels for MG hierarchy","HYPRE_StructPFMGSetMaxLevels",ex->max_levels,&ex->max_levels,PETSC_NULL);CHKERRQ(ierr); 1144 ierr = HYPRE_StructPFMGSetMaxLevels(ex->hsolver,ex->max_levels);CHKERRQ(ierr); 1145 1146 ierr = PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr); 1147 ierr = HYPRE_StructPFMGSetTol(ex->hsolver,ex->tol);CHKERRQ(ierr); 1148 ierr = PetscOptionsEList("-pc_pfmg_relax_type","Relax type for the up and down cycles","HYPRE_StructPFMGSetRelaxType",PFMGRelaxType,4,PFMGRelaxType[ex->relax_type],&ex->relax_type,PETSC_NULL);CHKERRQ(ierr); 1149 ierr = HYPRE_StructPFMGSetRelaxType(ex->hsolver, ex->relax_type);CHKERRQ(ierr); 1150 ierr = PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,2,PFMGRAPType[ex->rap_type],&ex->rap_type,PETSC_NULL);CHKERRQ(ierr); 1151 ierr = HYPRE_StructPFMGSetRAPType(ex->hsolver, ex->rap_type);CHKERRQ(ierr); 1152 ierr = PetscOptionsTail();CHKERRQ(ierr); 1153 PetscFunctionReturn(0); 1154 } 1155 1156 #undef __FUNCT__ 1157 #define __FUNCT__ "PCApply_PFMG" 1158 PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y) 1159 { 1160 PetscErrorCode ierr; 1161 PC_PFMG *ex = (PC_PFMG*) pc->data; 1162 PetscScalar *xx,*yy; 1163 int ilower[3],iupper[3]; 1164 Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data); 1165 1166 PetscFunctionBegin; 1167 ierr = DAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr); 1168 iupper[0] += ilower[0] - 1; 1169 iupper[1] += ilower[1] - 1; 1170 iupper[2] += ilower[2] - 1; 1171 1172 /* copy x values over to hypre */ 1173 ierr = HYPRE_StructVectorSetConstantValues(mx->hb,0.0);CHKERRQ(ierr); 1174 ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 1175 ierr = HYPRE_StructVectorSetBoxValues(mx->hb,ilower,iupper,xx);CHKERRQ(ierr); 1176 ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 1177 ierr = HYPRE_StructVectorAssemble(mx->hb);CHKERRQ(ierr); 1178 1179 ierr = HYPRE_StructPFMGSolve(ex->hsolver,mx->hmat,mx->hb,mx->hx);CHKERRQ(ierr); 1180 1181 /* copy solution values back to PETSc */ 1182 ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 1183 ierr = HYPRE_StructVectorGetBoxValues(mx->hx,ilower,iupper,yy);CHKERRQ(ierr); 1184 ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 1185 PetscFunctionReturn(0); 1186 } 1187 1188 #undef __FUNCT__ 1189 #define __FUNCT__ "PCApplyRichardson_PFMG" 1190 static PetscErrorCode PCApplyRichardson_PFMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason) 1191 { 1192 PC_PFMG *jac = (PC_PFMG*)pc->data; 1193 PetscErrorCode ierr; 1194 int oits; 1195 1196 PetscFunctionBegin; 1197 ierr = HYPRE_StructPFMGSetMaxIter(jac->hsolver,its*jac->its);CHKERRQ(ierr); 1198 ierr = HYPRE_StructPFMGSetTol(jac->hsolver,rtol);CHKERRQ(ierr); 1199 1200 ierr = PCApply_PFMG(pc,b,y);CHKERRQ(ierr); 1201 ierr = HYPRE_StructPFMGGetNumIterations(jac->hsolver,&oits);CHKERRQ(ierr); 1202 *outits = oits; 1203 if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS; 1204 else *reason = PCRICHARDSON_CONVERGED_RTOL; 1205 ierr = HYPRE_StructPFMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 1206 ierr = HYPRE_StructPFMGSetMaxIter(jac->hsolver,jac->its);CHKERRQ(ierr); 1207 PetscFunctionReturn(0); 1208 } 1209 1210 1211 #undef __FUNCT__ 1212 #define __FUNCT__ "PCSetUp_PFMG" 1213 PetscErrorCode PCSetUp_PFMG(PC pc) 1214 { 1215 PetscErrorCode ierr; 1216 PC_PFMG *ex = (PC_PFMG*) pc->data; 1217 Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data); 1218 PetscTruth flg; 1219 1220 PetscFunctionBegin; 1221 ierr = PetscTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg);CHKERRQ(ierr); 1222 if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner"); 1223 1224 /* create the hypre solver object and set its information */ 1225 if (ex->hsolver) { 1226 ierr = HYPRE_StructPFMGDestroy(ex->hsolver);CHKERRQ(ierr); 1227 } 1228 ierr = HYPRE_StructPFMGCreate(ex->hcomm,&ex->hsolver);CHKERRQ(ierr); 1229 ierr = PCSetFromOptions_PFMG(pc);CHKERRQ(ierr); 1230 ierr = HYPRE_StructPFMGSetup(ex->hsolver,mx->hmat,mx->hb,mx->hx);CHKERRQ(ierr); 1231 ierr = HYPRE_StructPFMGSetZeroGuess(ex->hsolver);CHKERRQ(ierr); 1232 1233 PetscFunctionReturn(0); 1234 } 1235 1236 1237 /*MC 1238 PCPFMG - the hypre PFMG multigrid solver 1239 1240 Level: advanced 1241 1242 Options Database: 1243 + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner 1244 . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid 1245 . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid 1246 . -pc_pfmg_tol <tol> tolerance of PFMG 1247 . -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 1248 - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin 1249 1250 Notes: This is for CELL-centered descretizations 1251 1252 This must be used with the MATHYPRESTRUCT matrix type. 1253 This is less general than in hypre, it supports only one block per process defined by a PETSc DA. 1254 1255 .seealso: PCMG, MATHYPRESTRUCT 1256 M*/ 1257 1258 EXTERN_C_BEGIN 1259 #undef __FUNCT__ 1260 #define __FUNCT__ "PCCreate_PFMG" 1261 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_PFMG(PC pc) 1262 { 1263 PetscErrorCode ierr; 1264 PC_PFMG *ex; 1265 1266 PetscFunctionBegin; 1267 ierr = PetscNew(PC_PFMG,&ex);CHKERRQ(ierr);\ 1268 pc->data = ex; 1269 1270 ex->its = 1; 1271 ex->tol = 1.e-8; 1272 ex->relax_type = 1; 1273 ex->rap_type = 0; 1274 ex->num_pre_relax = 1; 1275 ex->num_post_relax = 1; 1276 ex->max_levels = 0; 1277 1278 pc->ops->setfromoptions = PCSetFromOptions_PFMG; 1279 pc->ops->view = PCView_PFMG; 1280 pc->ops->destroy = PCDestroy_PFMG; 1281 pc->ops->apply = PCApply_PFMG; 1282 pc->ops->applyrichardson = PCApplyRichardson_PFMG; 1283 pc->ops->setup = PCSetUp_PFMG; 1284 ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr); 1285 ierr = HYPRE_StructPFMGCreate(ex->hcomm,&ex->hsolver);CHKERRQ(ierr); 1286 PetscFunctionReturn(0); 1287 } 1288 EXTERN_C_END 1289 1290 /* we know we are working with a HYPRE_SStructMatrix */ 1291 typedef struct { 1292 MPI_Comm hcomm; /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */ 1293 HYPRE_SStructSolver ss_solver; 1294 1295 /* keep copy of SYSPFMG options used so may view them */ 1296 int its; 1297 double tol; 1298 int relax_type; 1299 int num_pre_relax,num_post_relax; 1300 } PC_SysPFMG; 1301 1302 #undef __FUNCT__ 1303 #define __FUNCT__ "PCDestroy_SysPFMG" 1304 PetscErrorCode PCDestroy_SysPFMG(PC pc) 1305 { 1306 PetscErrorCode ierr; 1307 PC_SysPFMG *ex = (PC_SysPFMG*) pc->data; 1308 1309 PetscFunctionBegin; 1310 if (ex->ss_solver) {ierr = HYPRE_SStructSysPFMGDestroy(ex->ss_solver);CHKERRQ(ierr);} 1311 ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr); 1312 ierr = PetscFree(ex);CHKERRQ(ierr); 1313 PetscFunctionReturn(0); 1314 } 1315 1316 static const char *SysPFMGRelaxType[] = {"Weighted-Jacobi","Red/Black-Gauss-Seidel"}; 1317 1318 #undef __FUNCT__ 1319 #define __FUNCT__ "PCView_SysPFMG" 1320 PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer) 1321 { 1322 PetscErrorCode ierr; 1323 PetscTruth iascii; 1324 PC_SysPFMG *ex = (PC_SysPFMG*) pc->data; 1325 1326 PetscFunctionBegin; 1327 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1328 if (iascii) { 1329 ierr = PetscViewerASCIIPrintf(viewer," HYPRE SysPFMG preconditioning\n");CHKERRQ(ierr); 1330 ierr = PetscViewerASCIIPrintf(viewer," HYPRE SysPFMG: max iterations %d\n",ex->its);CHKERRQ(ierr); 1331 ierr = PetscViewerASCIIPrintf(viewer," HYPRE SysPFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr); 1332 ierr = PetscViewerASCIIPrintf(viewer," HYPRE SysPFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr); 1333 ierr = PetscViewerASCIIPrintf(viewer," HYPRE SysPFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr); 1334 } 1335 PetscFunctionReturn(0); 1336 } 1337 1338 1339 #undef __FUNCT__ 1340 #define __FUNCT__ "PCSetFromOptions_SysPFMG" 1341 PetscErrorCode PCSetFromOptions_SysPFMG(PC pc) 1342 { 1343 PetscErrorCode ierr; 1344 PC_SysPFMG *ex = (PC_SysPFMG*) pc->data; 1345 PetscTruth flg = PETSC_FALSE; 1346 1347 PetscFunctionBegin; 1348 ierr = PetscOptionsHead("SysPFMG options");CHKERRQ(ierr); 1349 ierr = PetscOptionsTruth("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 1350 if (flg) { 1351 int level=3; 1352 ierr = HYPRE_SStructSysPFMGSetPrintLevel(ex->ss_solver,level);CHKERRQ(ierr); 1353 } 1354 ierr = PetscOptionsInt("-pc_syspfmg_its","Number of iterations of SysPFMG to use as preconditioner","HYPRE_SStructSysPFMGSetMaxIter",ex->its,&ex->its,PETSC_NULL);CHKERRQ(ierr); 1355 ierr = HYPRE_SStructSysPFMGSetMaxIter(ex->ss_solver,ex->its);CHKERRQ(ierr); 1356 ierr = PetscOptionsInt("-pc_syspfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_SStructSysPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,PETSC_NULL);CHKERRQ(ierr); 1357 ierr = HYPRE_SStructSysPFMGSetNumPreRelax(ex->ss_solver,ex->num_pre_relax);CHKERRQ(ierr); 1358 ierr = PetscOptionsInt("-pc_syspfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_SStructSysPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,PETSC_NULL);CHKERRQ(ierr); 1359 ierr = HYPRE_SStructSysPFMGSetNumPostRelax(ex->ss_solver,ex->num_post_relax);CHKERRQ(ierr); 1360 1361 ierr = PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr); 1362 ierr = HYPRE_SStructSysPFMGSetTol(ex->ss_solver,ex->tol);CHKERRQ(ierr); 1363 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,PETSC_NULL);CHKERRQ(ierr); 1364 ierr = HYPRE_SStructSysPFMGSetRelaxType(ex->ss_solver, ex->relax_type);CHKERRQ(ierr); 1365 ierr = PetscOptionsTail();CHKERRQ(ierr); 1366 PetscFunctionReturn(0); 1367 } 1368 1369 #undef __FUNCT__ 1370 #define __FUNCT__ "PCApply_SysPFMG" 1371 PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y) 1372 { 1373 PetscErrorCode ierr; 1374 PC_SysPFMG *ex = (PC_SysPFMG*) pc->data; 1375 PetscScalar *xx,*yy; 1376 int ilower[3],iupper[3]; 1377 Mat_HYPRESStruct *mx = (Mat_HYPRESStruct *)(pc->pmat->data); 1378 int ordering= mx->dofs_order; 1379 int nvars= mx->nvars; 1380 int part= 0; 1381 int size; 1382 int i; 1383 1384 PetscFunctionBegin; 1385 ierr = DAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr); 1386 iupper[0] += ilower[0] - 1; 1387 iupper[1] += ilower[1] - 1; 1388 iupper[2] += ilower[2] - 1; 1389 1390 size= 1; 1391 for (i= 0; i< 3; i++) { 1392 size*= (iupper[i]-ilower[i]+1); 1393 } 1394 /* copy x values over to hypre for variable ordering */ 1395 if (ordering) { 1396 ierr = HYPRE_SStructVectorSetConstantValues(mx->ss_b,0.0);CHKERRQ(ierr); 1397 ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 1398 for (i= 0; i< nvars; i++) { 1399 ierr = HYPRE_SStructVectorSetBoxValues(mx->ss_b,part,ilower,iupper,i,xx+(size*i));CHKERRQ(ierr); 1400 } 1401 ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 1402 ierr = HYPRE_SStructVectorAssemble(mx->ss_b);CHKERRQ(ierr); 1403 1404 ierr = HYPRE_SStructMatrixMatvec(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x);CHKERRQ(ierr); 1405 ierr = HYPRE_SStructSysPFMGSolve(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x);CHKERRQ(ierr); 1406 1407 /* copy solution values back to PETSc */ 1408 ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 1409 for (i= 0; i< nvars; i++) { 1410 ierr = HYPRE_SStructVectorGetBoxValues(mx->ss_x,part,ilower,iupper,i,yy+(size*i));CHKERRQ(ierr); 1411 } 1412 ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 1413 } 1414 1415 else { /* nodal ordering must be mapped to variable ordering for sys_pfmg */ 1416 PetscScalar *z; 1417 int j, k; 1418 1419 ierr = PetscMalloc(nvars*size*sizeof(PetscScalar),&z);CHKERRQ(ierr); 1420 ierr = HYPRE_SStructVectorSetConstantValues(mx->ss_b,0.0);CHKERRQ(ierr); 1421 ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 1422 1423 /* transform nodal to hypre's variable ordering for sys_pfmg */ 1424 for (i= 0; i< size; i++) { 1425 k= i*nvars; 1426 for (j= 0; j< nvars; j++) { 1427 z[j*size+i]= xx[k+j]; 1428 } 1429 } 1430 for (i= 0; i< nvars; i++) { 1431 ierr = HYPRE_SStructVectorSetBoxValues(mx->ss_b,part,ilower,iupper,i,z+(size*i));CHKERRQ(ierr); 1432 } 1433 ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 1434 1435 ierr = HYPRE_SStructVectorAssemble(mx->ss_b);CHKERRQ(ierr); 1436 1437 ierr = HYPRE_SStructSysPFMGSolve(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x);CHKERRQ(ierr); 1438 1439 /* copy solution values back to PETSc */ 1440 ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 1441 for (i= 0; i< nvars; i++) { 1442 ierr = HYPRE_SStructVectorGetBoxValues(mx->ss_x,part,ilower,iupper,i,z+(size*i));CHKERRQ(ierr); 1443 } 1444 /* transform hypre's variable ordering for sys_pfmg to nodal ordering */ 1445 for (i= 0; i< size; i++) { 1446 k= i*nvars; 1447 for (j= 0; j< nvars; j++) { 1448 yy[k+j]= z[j*size+i]; 1449 } 1450 } 1451 ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 1452 1453 ierr = PetscFree(z);CHKERRQ(ierr); 1454 } 1455 1456 PetscFunctionReturn(0); 1457 } 1458 1459 #undef __FUNCT__ 1460 #define __FUNCT__ "PCApplyRichardson_SysPFMG" 1461 static PetscErrorCode PCApplyRichardson_SysPFMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason) 1462 { 1463 PC_SysPFMG *jac = (PC_SysPFMG*)pc->data; 1464 PetscErrorCode ierr; 1465 int oits; 1466 1467 PetscFunctionBegin; 1468 1469 ierr = HYPRE_SStructSysPFMGSetMaxIter(jac->ss_solver,its*jac->its);CHKERRQ(ierr); 1470 ierr = HYPRE_SStructSysPFMGSetTol(jac->ss_solver,rtol);CHKERRQ(ierr); 1471 1472 ierr = PCApply_SysPFMG(pc,b,y);CHKERRQ(ierr); 1473 ierr = HYPRE_SStructSysPFMGGetNumIterations(jac->ss_solver,&oits);CHKERRQ(ierr); 1474 *outits = oits; 1475 if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS; 1476 else *reason = PCRICHARDSON_CONVERGED_RTOL; 1477 ierr = HYPRE_SStructSysPFMGSetTol(jac->ss_solver,jac->tol);CHKERRQ(ierr); 1478 ierr = HYPRE_SStructSysPFMGSetMaxIter(jac->ss_solver,jac->its);CHKERRQ(ierr); 1479 PetscFunctionReturn(0); 1480 } 1481 1482 1483 #undef __FUNCT__ 1484 #define __FUNCT__ "PCSetUp_SysPFMG" 1485 PetscErrorCode PCSetUp_SysPFMG(PC pc) 1486 { 1487 PetscErrorCode ierr; 1488 PC_SysPFMG *ex = (PC_SysPFMG*) pc->data; 1489 Mat_HYPRESStruct *mx = (Mat_HYPRESStruct *)(pc->pmat->data); 1490 PetscTruth flg; 1491 1492 PetscFunctionBegin; 1493 ierr = PetscTypeCompare((PetscObject)pc->pmat,MATHYPRESSTRUCT,&flg);CHKERRQ(ierr); 1494 if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESSTRUCT with this preconditioner"); 1495 1496 /* create the hypre sstruct solver object and set its information */ 1497 if (ex->ss_solver) { 1498 ierr = HYPRE_SStructSysPFMGDestroy(ex->ss_solver);CHKERRQ(ierr); 1499 } 1500 ierr = HYPRE_SStructSysPFMGCreate(ex->hcomm,&ex->ss_solver);CHKERRQ(ierr); 1501 ierr = PCSetFromOptions_SysPFMG(pc);CHKERRQ(ierr); 1502 ierr = HYPRE_SStructSysPFMGSetZeroGuess(ex->ss_solver);CHKERRQ(ierr); 1503 ierr = HYPRE_SStructSysPFMGSetup(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x);CHKERRQ(ierr); 1504 1505 PetscFunctionReturn(0); 1506 } 1507 1508 1509 /*MC 1510 PCSysPFMG - the hypre SysPFMG multigrid solver 1511 1512 Level: advanced 1513 1514 Options Database: 1515 + -pc_syspfmg_its <its> number of iterations of SysPFMG to use as preconditioner 1516 . -pc_syspfmg_num_pre_relax <steps> number of smoothing steps before coarse grid 1517 . -pc_syspfmg_num_post_relax <steps> number of smoothing steps after coarse grid 1518 . -pc_syspfmg_tol <tol> tolerance of SysPFMG 1519 . -pc_syspfmg_relax_type -relaxation type for the up and down cycles, one of Weighted-Jacobi,Red/Black-Gauss-Seidel 1520 1521 Notes: This is for CELL-centered descretizations 1522 1523 This must be used with the MATHYPRESSTRUCT matrix type. 1524 This is less general than in hypre, it supports only one part, and one block per process defined by a PETSc DA. 1525 Also, only cell-centered variables. 1526 1527 .seealso: PCMG, MATHYPRESSTRUCT 1528 M*/ 1529 1530 EXTERN_C_BEGIN 1531 #undef __FUNCT__ 1532 #define __FUNCT__ "PCCreate_SysPFMG" 1533 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_SysPFMG(PC pc) 1534 { 1535 PetscErrorCode ierr; 1536 PC_SysPFMG *ex; 1537 1538 PetscFunctionBegin; 1539 ierr = PetscNew(PC_SysPFMG,&ex);CHKERRQ(ierr);\ 1540 pc->data = ex; 1541 1542 ex->its = 1; 1543 ex->tol = 1.e-8; 1544 ex->relax_type = 1; 1545 ex->num_pre_relax = 1; 1546 ex->num_post_relax = 1; 1547 1548 pc->ops->setfromoptions = PCSetFromOptions_SysPFMG; 1549 pc->ops->view = PCView_SysPFMG; 1550 pc->ops->destroy = PCDestroy_SysPFMG; 1551 pc->ops->apply = PCApply_SysPFMG; 1552 pc->ops->applyrichardson = PCApplyRichardson_SysPFMG; 1553 pc->ops->setup = PCSetUp_SysPFMG; 1554 ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr); 1555 ierr = HYPRE_SStructSysPFMGCreate(ex->hcomm,&ex->ss_solver);CHKERRQ(ierr); 1556 PetscFunctionReturn(0); 1557 } 1558 EXTERN_C_END 1559