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