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 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 139 PETSC_EXTERN PetscErrorCode PCRegister(const char[],PetscErrorCode(*)(PC)); 140 141 PETSC_EXTERN PetscErrorCode PCReset(PC); 142 PETSC_EXTERN PetscErrorCode PCDestroy(PC*); 143 PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC); 144 145 PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*); 146 PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*); 147 PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*); 148 149 PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat); 150 PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*); 151 PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *); 152 153 PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer); 154 PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer); 155 PETSC_STATIC_INLINE PetscErrorCode PCViewFromOptions(PC A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 156 157 PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]); 158 PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]); 159 PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]); 160 161 PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*); 162 163 /* 164 These are used to provide extra scaling of preconditioned 165 operator for time-stepping schemes like in SUNDIALS 166 */ 167 PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *); 168 PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec); 169 PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec); 170 PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec); 171 172 /* ------------- options specific to particular preconditioners --------- */ 173 174 /*E 175 PCJacobiType - What elements are used to form the Jacobi preconditioner 176 177 Level: intermediate 178 179 .seealso: 180 E*/ 181 typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType; 182 PETSC_EXTERN const char *const PCJacobiTypes[]; 183 184 PETSC_EXTERN PetscErrorCode PCJacobiSetType(PC,PCJacobiType); 185 PETSC_EXTERN PetscErrorCode PCJacobiGetType(PC,PCJacobiType*); 186 PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC,PetscBool); 187 PETSC_EXTERN PetscErrorCode PCJacobiGetUseAbs(PC,PetscBool*); 188 PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType); 189 PETSC_EXTERN PetscErrorCode PCSORGetSymmetric(PC,MatSORType*); 190 PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal); 191 PETSC_EXTERN PetscErrorCode PCSORGetOmega(PC,PetscReal*); 192 PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt); 193 PETSC_EXTERN PetscErrorCode PCSORGetIterations(PC,PetscInt*,PetscInt*); 194 195 PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal); 196 PETSC_EXTERN PetscErrorCode PCEisenstatGetOmega(PC,PetscReal*); 197 PETSC_EXTERN PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC,PetscBool); 198 PETSC_EXTERN PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC,PetscBool*); 199 200 PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]); 201 PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]); 202 PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]); 203 PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]); 204 205 PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec)); 206 PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)); 207 PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec)); 208 PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC)); 209 PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*)); 210 PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer)); 211 PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC)); 212 PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*); 213 PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**); 214 PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]); 215 PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]); 216 217 PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal); 218 219 PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType); 220 PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal); 221 222 PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage); 223 PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*); 224 PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC); 225 226 PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal); 227 PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal); 228 PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal); 229 230 PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType); 231 PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool ); 232 PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool ); 233 PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC,PetscBool); 234 PETSC_EXTERN PetscErrorCode PCFactorGetUseInPlace(PC,PetscBool*); 235 PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC,PetscBool); 236 PETSC_EXTERN PetscErrorCode PCFactorGetAllowDiagonalFill(PC,PetscBool*); 237 PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool); 238 239 PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt); 240 PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC,PetscInt*); 241 PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt); 242 243 PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]); 244 PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]); 245 PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt); 246 PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool); 247 PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*); 248 PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool); 249 250 /*E 251 PCASMType - Type of additive Schwarz method to use 252 253 $ PC_ASM_BASIC - Symmetric version where residuals from the ghost points are used 254 $ and computed values in ghost regions are added together. 255 $ Classical standard additive Schwarz. 256 $ PC_ASM_RESTRICT - Residuals from ghost points are used but computed values in ghost 257 $ region are discarded. 258 $ Default. 259 $ PC_ASM_INTERPOLATE - Residuals from ghost points are not used, computed values in ghost 260 $ region are added back in. 261 $ PC_ASM_NONE - Residuals from ghost points are not used, computed ghost values are 262 $ discarded. 263 $ Not very good. 264 265 Level: beginner 266 267 .seealso: PCASMSetType() 268 E*/ 269 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 270 PETSC_EXTERN const char *const PCASMTypes[]; 271 272 PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType); 273 PETSC_EXTERN PetscErrorCode PCASMGetType(PC,PCASMType*); 274 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]); 275 PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]); 276 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 277 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]); 278 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]); 279 280 /*E 281 PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain). 282 283 Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 284 domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 285 a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 286 (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 287 288 $ PC_GASM_BASIC - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 289 $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 290 $ from neighboring subdomains. 291 $ Classical standard additive Schwarz. 292 $ PC_GASM_RESTRICT - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 293 $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 294 $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 295 $ assumption). 296 $ Default. 297 $ PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 298 $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 299 $ from neighboring subdomains. 300 $ 301 $ PC_GASM_NONE - Residuals and corrections are zeroed out outside the local subdomains. 302 $ Not very good. 303 304 Level: beginner 305 306 .seealso: PCGASMSetType() 307 E*/ 308 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType; 309 PETSC_EXTERN const char *const PCGASMTypes[]; 310 311 PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]); 312 PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool); 313 PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt); 314 PETSC_EXTERN PetscErrorCode PCGASMSetDMSubdomains(PC,PetscBool); 315 PETSC_EXTERN PetscErrorCode PCGASMGetDMSubdomains(PC,PetscBool*); 316 PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool ); 317 318 PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType); 319 PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]); 320 PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]); 321 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 322 PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]); 323 PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]); 324 325 /*E 326 PCCompositeType - Determines how two or more preconditioner are composed 327 328 $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 329 $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 330 $ computed after the previous preconditioner application 331 $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 332 $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions) 333 $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 334 $ where first preconditioner is built from alpha I + S and second from 335 $ alpha I + R 336 337 Level: beginner 338 339 .seealso: PCCompositeSetType() 340 E*/ 341 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType; 342 PETSC_EXTERN const char *const PCCompositeTypes[]; 343 344 PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType); 345 PETSC_EXTERN PetscErrorCode PCCompositeGetType(PC,PCCompositeType*); 346 PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType); 347 PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *); 348 PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar); 349 350 PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt); 351 PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter); 352 PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*); 353 354 PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double); 355 PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt); 356 PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt); 357 PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt); 358 PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt); 359 PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt); 360 PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt); 361 PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt); 362 363 PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]); 364 PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]); 365 366 PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*); 367 PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType); 368 PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*); 369 PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt); 370 PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS); 371 PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*); 372 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool); 373 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*); 374 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC,PetscBool); 375 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC,PetscBool*); 376 PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC,PetscBool); 377 PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC,PetscBool*); 378 379 380 /*E 381 PCFieldSplitSchurPreType - Determines how to precondition Schur complement 382 383 Level: intermediate 384 385 .seealso: PCFieldSplitSetSchurPre() 386 E*/ 387 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; 388 PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[]; 389 390 /*E 391 PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 392 393 Level: intermediate 394 395 .seealso: PCFieldSplitSetSchurFactType() 396 E*/ 397 typedef enum { 398 PC_FIELDSPLIT_SCHUR_FACT_DIAG, 399 PC_FIELDSPLIT_SCHUR_FACT_LOWER, 400 PC_FIELDSPLIT_SCHUR_FACT_UPPER, 401 PC_FIELDSPLIT_SCHUR_FACT_FULL 402 } PCFieldSplitSchurFactType; 403 PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[]; 404 405 PETSC_EXTERN PETSC_DEPRECATED("Use PCFieldSplitSetSchurPre") PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat); 406 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC,PCFieldSplitSchurPreType,Mat); 407 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC,PCFieldSplitSchurPreType*,Mat*); 408 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType); 409 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*); 410 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC,Mat *S); 411 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC,Mat *S); 412 413 PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat); 414 PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat); 415 416 PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*); 417 418 PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]); 419 420 PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM); 421 PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*); 422 423 PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*); 424 PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*); 425 426 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal); 427 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt); 428 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool); 429 430 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal); 431 PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool); 432 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt); 433 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt); 434 /*E 435 PCPARMSGlobalType - Determines the global preconditioner method in PARMS 436 437 Level: intermediate 438 439 .seealso: PCPARMSSetGlobal() 440 E*/ 441 typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType; 442 PETSC_EXTERN const char *const PCPARMSGlobalTypes[]; 443 /*E 444 PCPARMSLocalType - Determines the local preconditioner method in PARMS 445 446 Level: intermediate 447 448 .seealso: PCPARMSSetLocal() 449 E*/ 450 typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType; 451 PETSC_EXTERN const char *const PCPARMSLocalTypes[]; 452 453 PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC,PCPARMSGlobalType); 454 PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC,PCPARMSLocalType); 455 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC,PetscReal,PetscInt); 456 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC,PetscInt); 457 PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC,PetscBool); 458 PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC,PetscInt,PetscInt,PetscInt); 459 460 /*E 461 PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method 462 463 Level: intermediate 464 465 .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseProl() 466 E*/ 467 typedef const char *PCGAMGType; 468 #define PCGAMGAGG "agg" 469 #define PCGAMGGEO "geo" 470 #define PCGAMGCLASSICAL "classical" 471 472 PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType); 473 PETSC_EXTERN PetscErrorCode PCGAMGGetType( PC,PCGAMGType*); 474 PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt); 475 PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool); 476 PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool); 477 PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt); 478 PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal); 479 PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt); 480 PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt); 481 PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC,PetscInt); 482 PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC,PetscBool); 483 PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool); 484 PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool); 485 PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void); 486 PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void); 487 488 typedef const char *PCGAMGClassicalType; 489 #define PCGAMGCLASSICALDIRECT "direct" 490 #define PCGAMGCLASSICALSTANDARD "standard" 491 PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC,PCGAMGClassicalType); 492 PETSC_EXTERN PetscErrorCode PCGAMGClassicalGetType(PC,PCGAMGClassicalType*); 493 494 #if defined(PETSC_HAVE_PCBDDC) 495 PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisMat(PC,Mat); 496 PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC,IS); 497 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt); 498 PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC,PetscInt); 499 PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace); 500 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS); 501 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC,IS); 502 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*); 503 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC,IS*); 504 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS); 505 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC,IS); 506 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*); 507 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC,IS*); 508 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]); 509 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC,PetscInt,IS[]); 510 PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode); 511 PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*); 512 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec); 513 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec); 514 #endif 515 516 PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool); 517 PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar); 518 PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec); 519 520 /*E 521 PCMGType - Determines the type of multigrid method that is run. 522 523 Level: beginner 524 525 Values: 526 + PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles() 527 . PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are 528 smoothed before updating the residual. This only uses the 529 down smoother, in the preconditioner the upper smoother is ignored 530 . PC_MG_FULL - same as multiplicative except one also performs grid sequencing, 531 that is starts on the coarsest grid, performs a cycle, interpolates 532 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 533 algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 534 calls the V-cycle only on the coarser level and has a post-smoother instead. 535 - PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level 536 from a finer 537 538 .seealso: PCMGSetType() 539 540 E*/ 541 typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType; 542 PETSC_EXTERN const char *const PCMGTypes[]; 543 #define PC_MG_CASCADE PC_MG_KASKADE; 544 545 /*E 546 PCMGCycleType - Use V-cycle or W-cycle 547 548 Level: beginner 549 550 Values: 551 + PC_MG_V_CYCLE 552 - PC_MG_W_CYCLE 553 554 .seealso: PCMGSetCycleType() 555 556 E*/ 557 typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType; 558 PETSC_EXTERN const char *const PCMGCycleTypes[]; 559 560 PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType); 561 PETSC_EXTERN PetscErrorCode PCMGGetType(PC,PCMGType*); 562 PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*); 563 PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*); 564 565 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt); 566 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt); 567 PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType); 568 PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType); 569 PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt); 570 PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt); 571 PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool); 572 PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*); 573 574 PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec); 575 PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec); 576 PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec); 577 578 PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat); 579 PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*); 580 PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat); 581 PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*); 582 PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec); 583 PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*); 584 PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat); 585 PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat,Vec,Vec,Vec); 586 587 /*E 588 PCExoticType - Face based or wirebasket based coarse grid space 589 590 Level: beginner 591 592 .seealso: PCExoticSetType(), PCEXOTIC 593 E*/ 594 typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType; 595 PETSC_EXTERN const char *const PCExoticTypes[]; 596 PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType); 597 598 #endif /* __PETSCPC_H */ 599