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