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