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