1 /* 2 Preconditioner module. 3 */ 4 #if !defined(__PETSCPC_H) 5 #define __PETSCPC_H 6 #include <petscmat.h> 7 #include <petscdmtypes.h> 8 9 PETSC_EXTERN PetscErrorCode PCInitializePackage(void); 10 11 /* 12 PCList contains the list of preconditioners currently registered 13 These are added with PCRegister() 14 */ 15 PETSC_EXTERN PetscFunctionList PCList; 16 17 /*S 18 PC - Abstract PETSc object that manages all preconditioners including direct solvers such as PCLU 19 20 Level: beginner 21 22 Concepts: preconditioners 23 24 .seealso: PCCreate(), PCSetType(), PCType (for list of available types) 25 S*/ 26 typedef struct _p_PC* PC; 27 28 /*J 29 PCType - String with the name of a PETSc preconditioner method. 30 31 Level: beginner 32 33 Notes: Click on the links below to see details on a particular solver 34 35 PCRegister() is used to register preconditioners that are then accessible via PCSetType() 36 37 .seealso: PCSetType(), PC, PCCreate(), PCRegister(), PCSetFromOptions() 38 J*/ 39 typedef const char* PCType; 40 #define PCNONE "none" 41 #define PCJACOBI "jacobi" 42 #define PCSOR "sor" 43 #define PCLU "lu" 44 #define PCSHELL "shell" 45 #define PCBJACOBI "bjacobi" 46 #define PCMG "mg" 47 #define PCEISENSTAT "eisenstat" 48 #define PCILU "ilu" 49 #define PCICC "icc" 50 #define PCASM "asm" 51 #define PCGASM "gasm" 52 #define PCKSP "ksp" 53 #define PCCOMPOSITE "composite" 54 #define PCREDUNDANT "redundant" 55 #define PCSPAI "spai" 56 #define PCNN "nn" 57 #define PCCHOLESKY "cholesky" 58 #define PCPBJACOBI "pbjacobi" 59 #define PCMAT "mat" 60 #define PCHYPRE "hypre" 61 #define PCPARMS "parms" 62 #define PCFIELDSPLIT "fieldsplit" 63 #define PCTFS "tfs" 64 #define PCML "ml" 65 #define PCGALERKIN "galerkin" 66 #define PCEXOTIC "exotic" 67 #define PCCP "cp" 68 #define PCBFBT "bfbt" 69 #define PCLSC "lsc" 70 #define PCPYTHON "python" 71 #define PCPFMG "pfmg" 72 #define PCSYSPFMG "syspfmg" 73 #define PCREDISTRIBUTE "redistribute" 74 #define PCSVD "svd" 75 #define PCGAMG "gamg" 76 #define PCSACUSP "sacusp" /* these four run on NVIDIA GPUs using CUSP */ 77 #define PCSACUSPPOLY "sacusppoly" 78 #define PCBICGSTABCUSP "bicgstabcusp" 79 #define PCAINVCUSP "ainvcusp" 80 #define PCBDDC "bddc" 81 #define PCKACZMARZ "kaczmarz" 82 83 /* Logging support */ 84 PETSC_EXTERN PetscClassId PC_CLASSID; 85 86 /*E 87 PCSide - If the preconditioner is to be applied to the left, right 88 or symmetrically around the operator. 89 90 Level: beginner 91 92 .seealso: 93 E*/ 94 typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide; 95 #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 96 PETSC_EXTERN const char *const *const PCSides; 97 98 PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*); 99 PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType); 100 PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*); 101 PETSC_EXTERN PetscErrorCode PCSetUp(PC); 102 PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC); 103 PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec); 104 PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec); 105 PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec); 106 PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec); 107 PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec); 108 PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *); 109 PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec); 110 PETSC_EXTERN PetscErrorCode PCSetReusePreconditioner(PC,PetscBool); 111 PETSC_EXTERN PetscErrorCode PCGetReusePreconditioner(PC,PetscBool*); 112 113 #define PC_FILE_CLASSID 1211222 114 115 /*E 116 PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates 117 118 Level: advanced 119 120 Notes: this must match petsc-finclude/petscpc.h and the KSPConvergedReason values in petscksp.h 121 122 .seealso: PCApplyRichardson() 123 E*/ 124 typedef enum { 125 PCRICHARDSON_CONVERGED_RTOL = 2, 126 PCRICHARDSON_CONVERGED_ATOL = 3, 127 PCRICHARDSON_CONVERGED_ITS = 4, 128 PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason; 129 130 PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*); 131 PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *); 132 PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool); 133 PETSC_EXTERN PetscErrorCode PCGetInitialGuessNonzero(PC,PetscBool*); 134 PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool); 135 PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*); 136 137 PETSC_EXTERN PetscErrorCode PCRegisterAll(void); 138 PETSC_EXTERN PetscBool PCRegisterAllCalled; 139 140 PETSC_EXTERN PetscErrorCode PCRegister(const char[],PetscErrorCode(*)(PC)); 141 142 PETSC_EXTERN PetscErrorCode PCReset(PC); 143 PETSC_EXTERN PetscErrorCode PCDestroy(PC*); 144 PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC); 145 146 PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*); 147 PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*); 148 PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*); 149 150 PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat); 151 PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*); 152 PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *); 153 154 PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer); 155 PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer); 156 PETSC_STATIC_INLINE PetscErrorCode PCViewFromOptions(PC A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 157 158 PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]); 159 PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]); 160 PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]); 161 162 PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*); 163 164 /* 165 These are used to provide extra scaling of preconditioned 166 operator for time-stepping schemes like in SUNDIALS 167 */ 168 PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *); 169 PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec); 170 PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec); 171 PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec); 172 173 /* ------------- options specific to particular preconditioners --------- */ 174 175 /*E 176 PCJacobiType - What elements are used to form the Jacobi preconditioner 177 178 Level: intermediate 179 180 .seealso: 181 E*/ 182 typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType; 183 PETSC_EXTERN const char *const PCJacobiTypes[]; 184 185 PETSC_EXTERN PetscErrorCode PCJacobiSetType(PC,PCJacobiType); 186 PETSC_EXTERN PetscErrorCode PCJacobiGetType(PC,PCJacobiType*); 187 PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC,PetscBool); 188 PETSC_EXTERN PetscErrorCode PCJacobiGetUseAbs(PC,PetscBool*); 189 PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType); 190 PETSC_EXTERN PetscErrorCode PCSORGetSymmetric(PC,MatSORType*); 191 PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal); 192 PETSC_EXTERN PetscErrorCode PCSORGetOmega(PC,PetscReal*); 193 PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt); 194 PETSC_EXTERN PetscErrorCode PCSORGetIterations(PC,PetscInt*,PetscInt*); 195 196 PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal); 197 PETSC_EXTERN PetscErrorCode PCEisenstatGetOmega(PC,PetscReal*); 198 PETSC_EXTERN PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC,PetscBool); 199 PETSC_EXTERN PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC,PetscBool*); 200 201 PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]); 202 PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]); 203 PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]); 204 PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]); 205 206 PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec)); 207 PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)); 208 PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec)); 209 PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC)); 210 PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*)); 211 PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer)); 212 PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC)); 213 PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*); 214 PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**); 215 PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]); 216 PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]); 217 218 PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal); 219 220 PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType); 221 PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal); 222 223 PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage); 224 PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*); 225 PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC); 226 227 PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal); 228 PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal); 229 PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal); 230 231 PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType); 232 PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool ); 233 PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool ); 234 PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC,PetscBool); 235 PETSC_EXTERN PetscErrorCode PCFactorGetUseInPlace(PC,PetscBool*); 236 PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC,PetscBool); 237 PETSC_EXTERN PetscErrorCode PCFactorGetAllowDiagonalFill(PC,PetscBool*); 238 PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool); 239 240 PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt); 241 PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC,PetscInt*); 242 PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt); 243 244 PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]); 245 PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]); 246 PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt); 247 PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool); 248 PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*); 249 PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool); 250 251 /*E 252 PCASMType - Type of additive Schwarz method to use 253 254 $ PC_ASM_BASIC - Symmetric version where residuals from the ghost points are used 255 $ and computed values in ghost regions are added together. 256 $ Classical standard additive Schwarz. 257 $ PC_ASM_RESTRICT - Residuals from ghost points are used but computed values in ghost 258 $ region are discarded. 259 $ Default. 260 $ PC_ASM_INTERPOLATE - Residuals from ghost points are not used, computed values in ghost 261 $ region are added back in. 262 $ PC_ASM_NONE - Residuals from ghost points are not used, computed ghost values are 263 $ discarded. 264 $ Not very good. 265 266 Level: beginner 267 268 .seealso: PCASMSetType() 269 E*/ 270 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 271 PETSC_EXTERN const char *const PCASMTypes[]; 272 273 PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType); 274 PETSC_EXTERN PetscErrorCode PCASMGetType(PC,PCASMType*); 275 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]); 276 PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]); 277 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 278 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]); 279 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]); 280 281 /*E 282 PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain). 283 284 Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 285 domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 286 a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 287 (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 288 289 $ PC_GASM_BASIC - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 290 $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 291 $ from neighboring subdomains. 292 $ Classical standard additive Schwarz. 293 $ PC_GASM_RESTRICT - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 294 $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 295 $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 296 $ assumption). 297 $ Default. 298 $ PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 299 $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 300 $ from neighboring subdomains. 301 $ 302 $ PC_GASM_NONE - Residuals and corrections are zeroed out outside the local subdomains. 303 $ Not very good. 304 305 Level: beginner 306 307 .seealso: PCGASMSetType() 308 E*/ 309 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType; 310 PETSC_EXTERN const char *const PCGASMTypes[]; 311 312 PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]); 313 PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool); 314 PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt); 315 PETSC_EXTERN PetscErrorCode PCGASMSetDMSubdomains(PC,PetscBool); 316 PETSC_EXTERN PetscErrorCode PCGASMGetDMSubdomains(PC,PetscBool*); 317 PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool ); 318 319 PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType); 320 PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]); 321 PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]); 322 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 323 PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]); 324 PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]); 325 326 /*E 327 PCCompositeType - Determines how two or more preconditioner are composed 328 329 $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 330 $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 331 $ computed after the previous preconditioner application 332 $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 333 $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions) 334 $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 335 $ where first preconditioner is built from alpha I + S and second from 336 $ alpha I + R 337 338 Level: beginner 339 340 .seealso: PCCompositeSetType() 341 E*/ 342 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType; 343 PETSC_EXTERN const char *const PCCompositeTypes[]; 344 345 PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType); 346 PETSC_EXTERN PetscErrorCode PCCompositeGetType(PC,PCCompositeType*); 347 PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType); 348 PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *); 349 PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar); 350 351 PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt); 352 PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter); 353 PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*); 354 355 PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double); 356 PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt); 357 PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt); 358 PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt); 359 PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt); 360 PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt); 361 PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt); 362 PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt); 363 364 PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]); 365 PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]); 366 PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteGradient(PC,Mat); 367 PETSC_EXTERN PetscErrorCode PCHYPRESetEdgeConstantVectors(PC,Vec,Vec,Vec); 368 PETSC_EXTERN PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC,Mat); 369 PETSC_EXTERN PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC,Mat); 370 PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]); 371 PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]); 372 373 PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*); 374 PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType); 375 PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*); 376 PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt); 377 PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS); 378 PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*); 379 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool); 380 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*); 381 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC,PetscBool); 382 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC,PetscBool*); 383 PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC,PetscBool); 384 PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC,PetscBool*); 385 386 387 /*E 388 PCFieldSplitSchurPreType - Determines how to precondition Schur complement 389 390 Level: intermediate 391 392 .seealso: PCFieldSplitSetSchurPre() 393 E*/ 394 typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_SELFP,PC_FIELDSPLIT_SCHUR_PRE_A11,PC_FIELDSPLIT_SCHUR_PRE_USER,PC_FIELDSPLIT_SCHUR_PRE_FULL} PCFieldSplitSchurPreType; 395 PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[]; 396 397 /*E 398 PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 399 400 Level: intermediate 401 402 .seealso: PCFieldSplitSetSchurFactType() 403 E*/ 404 typedef enum { 405 PC_FIELDSPLIT_SCHUR_FACT_DIAG, 406 PC_FIELDSPLIT_SCHUR_FACT_LOWER, 407 PC_FIELDSPLIT_SCHUR_FACT_UPPER, 408 PC_FIELDSPLIT_SCHUR_FACT_FULL 409 } PCFieldSplitSchurFactType; 410 PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[]; 411 412 PETSC_EXTERN PETSC_DEPRECATED("Use PCFieldSplitSetSchurPre") PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat); 413 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC,PCFieldSplitSchurPreType,Mat); 414 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC,PCFieldSplitSchurPreType*,Mat*); 415 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType); 416 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*); 417 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC,Mat *S); 418 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC,Mat *S); 419 420 PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat); 421 PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat); 422 423 PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*); 424 425 PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]); 426 427 PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM); 428 PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*); 429 430 PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*); 431 PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*); 432 433 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal); 434 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt); 435 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool); 436 437 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal); 438 PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool); 439 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt); 440 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt); 441 /*E 442 PCPARMSGlobalType - Determines the global preconditioner method in PARMS 443 444 Level: intermediate 445 446 .seealso: PCPARMSSetGlobal() 447 E*/ 448 typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType; 449 PETSC_EXTERN const char *const PCPARMSGlobalTypes[]; 450 /*E 451 PCPARMSLocalType - Determines the local preconditioner method in PARMS 452 453 Level: intermediate 454 455 .seealso: PCPARMSSetLocal() 456 E*/ 457 typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType; 458 PETSC_EXTERN const char *const PCPARMSLocalTypes[]; 459 460 PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC,PCPARMSGlobalType); 461 PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC,PCPARMSLocalType); 462 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC,PetscReal,PetscInt); 463 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC,PetscInt); 464 PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC,PetscBool); 465 PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC,PetscInt,PetscInt,PetscInt); 466 467 /*E 468 PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method 469 470 Level: intermediate 471 472 .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseInterpolation() 473 E*/ 474 typedef const char *PCGAMGType; 475 #define PCGAMGAGG "agg" 476 #define PCGAMGGEO "geo" 477 #define PCGAMGCLASSICAL "classical" 478 479 PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType); 480 PETSC_EXTERN PetscErrorCode PCGAMGGetType( PC,PCGAMGType*); 481 PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt); 482 PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool); 483 PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool); 484 PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt); 485 PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal); 486 PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt); 487 PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt); 488 PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC,PetscInt); 489 PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC,PetscBool); 490 PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool); 491 PETSC_EXTERN PetscErrorCode PCGAMGSetReuseInterpolation(PC,PetscBool); 492 PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void); 493 PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void); 494 PETSC_EXTERN PetscErrorCode PCGAMGRegister(PCGAMGType,PetscErrorCode (*)(PC)); 495 496 typedef const char *PCGAMGClassicalType; 497 #define PCGAMGCLASSICALDIRECT "direct" 498 #define PCGAMGCLASSICALSTANDARD "standard" 499 PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC,PCGAMGClassicalType); 500 PETSC_EXTERN PetscErrorCode PCGAMGClassicalGetType(PC,PCGAMGClassicalType*); 501 502 PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisMat(PC,Mat); 503 PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC,IS); 504 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt); 505 PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC,PetscInt); 506 PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace); 507 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS); 508 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC,IS); 509 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*); 510 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC,IS*); 511 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS); 512 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC,IS); 513 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*); 514 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC,IS*); 515 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]); 516 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC,PetscInt,IS[]); 517 PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode); 518 PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*); 519 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec); 520 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec); 521 522 PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool); 523 PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar); 524 PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec); 525 526 /*E 527 PCMGType - Determines the type of multigrid method that is run. 528 529 Level: beginner 530 531 Values: 532 + PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles() 533 . PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are 534 smoothed before updating the residual. This only uses the 535 down smoother, in the preconditioner the upper smoother is ignored 536 . PC_MG_FULL - same as multiplicative except one also performs grid sequencing, 537 that is starts on the coarsest grid, performs a cycle, interpolates 538 to the next, performs a cycle etc. This is much like the F-cycle presented in "Multigrid" by Trottenberg, Oosterlee, Schuller page 49, but that 539 algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 540 calls the V-cycle only on the coarser level and has a post-smoother instead. 541 - PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level 542 from a finer 543 544 .seealso: PCMGSetType() 545 546 E*/ 547 typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType; 548 PETSC_EXTERN const char *const PCMGTypes[]; 549 #define PC_MG_CASCADE PC_MG_KASKADE; 550 551 /*E 552 PCMGCycleType - Use V-cycle or W-cycle 553 554 Level: beginner 555 556 Values: 557 + PC_MG_V_CYCLE 558 - PC_MG_W_CYCLE 559 560 .seealso: PCMGSetCycleType() 561 562 E*/ 563 typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType; 564 PETSC_EXTERN const char *const PCMGCycleTypes[]; 565 566 PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType); 567 PETSC_EXTERN PetscErrorCode PCMGGetType(PC,PCMGType*); 568 PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*); 569 PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*); 570 571 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt); 572 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt); 573 PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType); 574 PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType); 575 PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt); 576 PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt); 577 PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool); 578 PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*); 579 580 PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec); 581 PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec); 582 PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec); 583 584 PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat); 585 PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*); 586 PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat); 587 PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*); 588 PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec); 589 PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*); 590 PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat); 591 PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat,Vec,Vec,Vec); 592 593 /*E 594 PCExoticType - Face based or wirebasket based coarse grid space 595 596 Level: beginner 597 598 .seealso: PCExoticSetType(), PCEXOTIC 599 E*/ 600 typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType; 601 PETSC_EXTERN const char *const PCExoticTypes[]; 602 PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType); 603 604 #endif /* __PETSCPC_H */ 605