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(const char[]); 10 11 /* 12 PCList contains the list of preconditioners currently registered 13 These are added with the PCRegisterDynamic() macro 14 */ 15 PETSC_EXTERN PetscFunctionList PCList; 16 17 /*S 18 PC - Abstract PETSc object that manages all preconditioners 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 or the creation function 30 with an optional dynamic library name, for example 31 http://www.mcs.anl.gov/petsc/lib.a:mypccreate() 32 33 Level: beginner 34 35 Notes: Click on the links below to see details on a particular solver 36 37 .seealso: PCSetType(), PC, PCCreate() 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 PCHMPI "hmpi" 68 #define PCSUPPORTGRAPH "supportgraph" 69 #define PCASA "asa" 70 #define PCCP "cp" 71 #define PCBFBT "bfbt" 72 #define PCLSC "lsc" 73 #define PCPYTHON "python" 74 #define PCPFMG "pfmg" 75 #define PCSYSPFMG "syspfmg" 76 #define PCREDISTRIBUTE "redistribute" 77 #define PCSVD "svd" 78 #define PCGAMG "gamg" 79 #define PCSACUSP "sacusp" /* these four run on NVIDIA GPUs using CUSP */ 80 #define PCSACUSPPOLY "sacusppoly" 81 #define PCBICGSTABCUSP "bicgstabcusp" 82 #define PCAINVCUSP "ainvcusp" 83 #define PCBDDC "bddc" 84 85 /* Logging support */ 86 PETSC_EXTERN PetscClassId PC_CLASSID; 87 88 /*E 89 PCSide - If the preconditioner is to be applied to the left, right 90 or symmetrically around the operator. 91 92 Level: beginner 93 94 .seealso: 95 E*/ 96 typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide; 97 #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 98 PETSC_EXTERN const char *const *const PCSides; 99 100 PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*); 101 PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType); 102 PETSC_EXTERN PetscErrorCode PCSetUp(PC); 103 PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC); 104 PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec); 105 PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec); 106 PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec); 107 PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec); 108 PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec); 109 PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *); 110 PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec); 111 112 #define PC_FILE_CLASSID 1211222 113 114 /*E 115 PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates 116 117 Level: advanced 118 119 Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h 120 121 .seealso: PCApplyRichardson() 122 E*/ 123 typedef enum { 124 PCRICHARDSON_CONVERGED_RTOL = 2, 125 PCRICHARDSON_CONVERGED_ATOL = 3, 126 PCRICHARDSON_CONVERGED_ITS = 4, 127 PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason; 128 129 PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*); 130 PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *); 131 PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool ); 132 133 PETSC_EXTERN PetscErrorCode PCRegisterDestroy(void); 134 PETSC_EXTERN PetscErrorCode PCRegisterAll(const char[]); 135 PETSC_EXTERN PetscBool PCRegisterAllCalled; 136 137 PETSC_EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC)); 138 139 /*MC 140 PCRegisterDynamic - Adds a method to the preconditioner package. 141 142 Synopsis: 143 #include "petscpc.h" 144 PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC)) 145 146 Not collective 147 148 Input Parameters: 149 + name_solver - name of a new user-defined solver 150 . path - path (either absolute or relative) the library containing this solver 151 . name_create - name of routine to create method context 152 - routine_create - routine to create method context 153 154 Notes: 155 PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners. 156 157 If dynamic libraries are used, then the fourth input argument (routine_create) 158 is ignored. 159 160 Sample usage: 161 .vb 162 PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib", 163 "MySolverCreate",MySolverCreate); 164 .ve 165 166 Then, your solver can be chosen with the procedural interface via 167 $ PCSetType(pc,"my_solver") 168 or at runtime via the option 169 $ -pc_type my_solver 170 171 Level: advanced 172 173 Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, or ${any environmental variable} 174 occuring in pathname will be replaced with appropriate values. 175 If your function is not being put into a shared library then use PCRegister() instead 176 177 .keywords: PC, register 178 179 .seealso: PCRegisterAll(), PCRegisterDestroy() 180 M*/ 181 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 182 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0) 183 #else 184 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d) 185 #endif 186 187 PETSC_EXTERN PetscErrorCode PCReset(PC); 188 PETSC_EXTERN PetscErrorCode PCDestroy(PC*); 189 PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC); 190 PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*); 191 192 PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*); 193 PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*); 194 PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*); 195 196 PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure); 197 PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*); 198 PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *); 199 200 PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer); 201 PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer); 202 203 PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]); 204 PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]); 205 PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]); 206 207 PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*); 208 209 /* 210 These are used to provide extra scaling of preconditioned 211 operator for time-stepping schemes like in SUNDIALS 212 */ 213 PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *); 214 PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec); 215 PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec); 216 PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec); 217 218 /* ------------- options specific to particular preconditioners --------- */ 219 220 PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC); 221 PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC); 222 PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC); 223 PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType); 224 PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal); 225 PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt); 226 227 PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal); 228 PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC); 229 230 #define USE_PRECONDITIONER_MATRIX 0 231 #define USE_TRUE_MATRIX 1 232 PETSC_EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC); 233 PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]); 234 PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]); 235 236 PETSC_EXTERN PetscErrorCode PCKSPSetUseTrue(PC); 237 238 PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec)); 239 PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)); 240 PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec)); 241 PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC)); 242 PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*)); 243 PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer)); 244 PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC)); 245 PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**); 246 PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*); 247 PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]); 248 PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]); 249 250 PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal); 251 252 PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType); 253 PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal); 254 255 PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage); 256 PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*); 257 PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC); 258 259 PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal); 260 PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal); 261 PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal); 262 263 PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType); 264 PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool ); 265 PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool ); 266 PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC); 267 PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC); 268 PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool ); 269 270 PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt); 271 PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt); 272 273 PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]); 274 PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]); 275 PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt); 276 PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool); 277 PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*); 278 PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool ); 279 280 /*E 281 PCASMType - Type of additive Schwarz method to use 282 283 $ PC_ASM_BASIC - Symmetric version where residuals from the ghost points are used 284 $ and computed values in ghost regions are added together. 285 $ Classical standard additive Schwarz. 286 $ PC_ASM_RESTRICT - Residuals from ghost points are used but computed values in ghost 287 $ region are discarded. 288 $ Default. 289 $ PC_ASM_INTERPOLATE - Residuals from ghost points are not used, computed values in ghost 290 $ region are added back in. 291 $ PC_ASM_NONE - Residuals from ghost points are not used, computed ghost values are 292 $ discarded. 293 $ Not very good. 294 295 Level: beginner 296 297 .seealso: PCASMSetType() 298 E*/ 299 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 300 PETSC_EXTERN const char *const PCASMTypes[]; 301 302 PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType); 303 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]); 304 PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]); 305 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 306 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]); 307 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]); 308 309 /*E 310 PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain). 311 312 Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 313 domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 314 a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 315 (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 316 317 $ PC_GASM_BASIC - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 318 $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 319 $ from neighboring subdomains. 320 $ Classical standard additive Schwarz. 321 $ PC_GASM_RESTRICT - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 322 $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 323 $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 324 $ assumption). 325 $ Default. 326 $ PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 327 $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 328 $ from neighboring subdomains. 329 $ 330 $ PC_GASM_NONE - Residuals and corrections are zeroed out outside the local subdomains. 331 $ Not very good. 332 333 Level: beginner 334 335 .seealso: PCGASMSetType() 336 E*/ 337 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType; 338 PETSC_EXTERN const char *const PCGASMTypes[]; 339 340 PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]); 341 PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool); 342 PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt); 343 PETSC_EXTERN PetscErrorCode PCGASMSetDMSubdomains(PC,PetscBool); 344 PETSC_EXTERN PetscErrorCode PCGASMGetDMSubdomains(PC,PetscBool*); 345 PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool ); 346 347 PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType); 348 PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]); 349 PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]); 350 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 351 PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]); 352 PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]); 353 354 /*E 355 PCCompositeType - Determines how two or more preconditioner are composed 356 357 $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 358 $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 359 $ computed after the previous preconditioner application 360 $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 361 $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions) 362 $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 363 $ where first preconditioner is built from alpha I + S and second from 364 $ alpha I + R 365 366 Level: beginner 367 368 .seealso: PCCompositeSetType() 369 E*/ 370 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType; 371 PETSC_EXTERN const char *const PCCompositeTypes[]; 372 373 PETSC_EXTERN PetscErrorCode PCCompositeSetUseTrue(PC); 374 PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType); 375 PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType); 376 PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *); 377 PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar); 378 379 PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt); 380 PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter); 381 PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*); 382 383 PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double); 384 PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt); 385 PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt); 386 PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt); 387 PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt); 388 PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt); 389 PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt); 390 PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt); 391 392 PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]); 393 PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]); 394 PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]); 395 PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]); 396 397 PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*); 398 PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*); 399 PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType); 400 PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt); 401 PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS); 402 PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*); 403 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool); 404 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*); 405 406 /*E 407 PCFieldSplitSchurPreType - Determines how to precondition Schur complement 408 409 Level: intermediate 410 411 .seealso: PCFieldSplitSchurPrecondition() 412 E*/ 413 typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_A11,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType; 414 PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[]; 415 416 /*E 417 PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 418 419 Level: intermediate 420 421 .seealso: PCFieldSplitSetSchurFactType() 422 E*/ 423 typedef enum { 424 PC_FIELDSPLIT_SCHUR_FACT_DIAG, 425 PC_FIELDSPLIT_SCHUR_FACT_LOWER, 426 PC_FIELDSPLIT_SCHUR_FACT_UPPER, 427 PC_FIELDSPLIT_SCHUR_FACT_FULL 428 } PCFieldSplitSchurFactType; 429 PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[]; 430 431 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat); 432 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType); 433 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*); 434 435 PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat); 436 PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat); 437 438 PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*); 439 PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *); 440 441 PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]); 442 443 PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM); 444 PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*); 445 446 PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*); 447 PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*); 448 449 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal); 450 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt); 451 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool); 452 453 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal); 454 PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool); 455 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt); 456 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt); 457 /*E 458 PCPARMSGlobalType - Determines the global preconditioner method in PARMS 459 460 Level: intermediate 461 462 .seealso: PCPARMSSetGlobal() 463 E*/ 464 typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType; 465 PETSC_EXTERN const char *const PCPARMSGlobalTypes[]; 466 /*E 467 PCPARMSLocalType - Determines the local preconditioner method in PARMS 468 469 Level: intermediate 470 471 .seealso: PCPARMSSetLocal() 472 E*/ 473 typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType; 474 PETSC_EXTERN const char *const PCPARMSLocalTypes[]; 475 476 PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type); 477 PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type); 478 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits); 479 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart); 480 PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym); 481 PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2); 482 483 PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt); 484 PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool); 485 PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool); 486 PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt); 487 PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal); 488 PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt); 489 PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt); 490 typedef const char* PCGAMGType; 491 PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType ); 492 PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n); 493 PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n); 494 PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool); 495 PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool); 496 497 #if defined(PETSC_HAVE_PCBDDC) 498 /* Enum defining how to treat the coarse problem */ 499 typedef enum {SEQUENTIAL_BDDC,REPLICATED_BDDC,PARALLEL_BDDC,MULTILEVEL_BDDC} CoarseProblemType; 500 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt); 501 PETSC_EXTERN PetscErrorCode PCBDDCSetMaxLevels(PC,PetscInt); 502 PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace); 503 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS); 504 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*); 505 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS); 506 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*); 507 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseProblemType(PC,CoarseProblemType); 508 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(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 PCMGSetLevels(PC,PetscInt,MPI_Comm*); 561 PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*); 562 563 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt); 564 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt); 565 PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType); 566 PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType); 567 PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt); 568 PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt); 569 PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool); 570 PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*); 571 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 PCMGDefaultResidual(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