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