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