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