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