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