1 /* 2 Preconditioner module. 3 */ 4 #if !defined(__PETSCPC_H) 5 #define __PETSCPC_H 6 #include "petscdm.h" 7 PETSC_EXTERN_CXX_BEGIN 8 9 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 extern PetscFList 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 /*E 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 E*/ 39 #define PCType char* 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 PCPROMETHEUS "prometheus" 66 #define PCGALERKIN "galerkin" 67 #define PCEXOTIC "exotic" 68 #define PCOPENMP "openmp" 69 #define PCSUPPORTGRAPH "supportgraph" 70 #define PCASA "asa" 71 #define PCCP "cp" 72 #define PCBFBT "bfbt" 73 #define PCLSC "lsc" 74 #define PCPYTHON "python" 75 #define PCPFMG "pfmg" 76 #define PCSYSPFMG "syspfmg" 77 #define PCREDISTRIBUTE "redistribute" 78 #define PCSACUSP "sacusp" 79 #define PCSACUSPPOLY "sacusppoly" 80 #define PCBICGSTABCUSP "bicgstabcusp" 81 #define PCSVD "svd" 82 #define PCAINVCUSP "ainvcusp" 83 84 /* Logging support */ 85 extern PetscClassId PC_CLASSID; 86 87 /*E 88 PCSide - If the preconditioner is to be applied to the left, right 89 or symmetrically around the operator. 90 91 Level: beginner 92 93 .seealso: 94 E*/ 95 typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide; 96 extern const char *PCSides[]; 97 98 extern PetscErrorCode PCCreate(MPI_Comm,PC*); 99 extern PetscErrorCode PCSetType(PC,const PCType); 100 extern PetscErrorCode PCSetUp(PC); 101 extern PetscErrorCode PCSetUpOnBlocks(PC); 102 extern PetscErrorCode PCApply(PC,Vec,Vec); 103 extern PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec); 104 extern PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec); 105 extern PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec); 106 extern PetscErrorCode PCApplyTranspose(PC,Vec,Vec); 107 extern PetscErrorCode PCApplyTransposeExists(PC,PetscBool *); 108 extern PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec); 109 110 /*E 111 PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates 112 113 Level: advanced 114 115 Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h 116 117 .seealso: PCApplyRichardson() 118 E*/ 119 typedef enum { 120 PCRICHARDSON_CONVERGED_RTOL = 2, 121 PCRICHARDSON_CONVERGED_ATOL = 3, 122 PCRICHARDSON_CONVERGED_ITS = 4, 123 PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason; 124 125 extern PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*); 126 extern PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *); 127 extern PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool ); 128 129 extern PetscErrorCode PCRegisterDestroy(void); 130 extern PetscErrorCode PCRegisterAll(const char[]); 131 extern PetscBool PCRegisterAllCalled; 132 133 extern PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC)); 134 135 /*MC 136 PCRegisterDynamic - Adds a method to the preconditioner package. 137 138 Synopsis: 139 PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC)) 140 141 Not collective 142 143 Input Parameters: 144 + name_solver - name of a new user-defined solver 145 . path - path (either absolute or relative) the library containing this solver 146 . name_create - name of routine to create method context 147 - routine_create - routine to create method context 148 149 Notes: 150 PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners. 151 152 If dynamic libraries are used, then the fourth input argument (routine_create) 153 is ignored. 154 155 Sample usage: 156 .vb 157 PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib", 158 "MySolverCreate",MySolverCreate); 159 .ve 160 161 Then, your solver can be chosen with the procedural interface via 162 $ PCSetType(pc,"my_solver") 163 or at runtime via the option 164 $ -pc_type my_solver 165 166 Level: advanced 167 168 Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, or ${any environmental variable} 169 occuring in pathname will be replaced with appropriate values. 170 If your function is not being put into a shared library then use PCRegister() instead 171 172 .keywords: PC, register 173 174 .seealso: PCRegisterAll(), PCRegisterDestroy() 175 M*/ 176 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 177 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0) 178 #else 179 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d) 180 #endif 181 182 extern PetscErrorCode PCReset(PC); 183 extern PetscErrorCode PCDestroy(PC*); 184 extern PetscErrorCode PCSetFromOptions(PC); 185 extern PetscErrorCode PCGetType(PC,const PCType*); 186 187 extern PetscErrorCode PCFactorGetMatrix(PC,Mat*); 188 extern PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*); 189 extern PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*); 190 191 extern PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure); 192 extern PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*); 193 extern PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *); 194 195 extern PetscErrorCode PCView(PC,PetscViewer); 196 197 extern PetscErrorCode PCSetOptionsPrefix(PC,const char[]); 198 extern PetscErrorCode PCAppendOptionsPrefix(PC,const char[]); 199 extern PetscErrorCode PCGetOptionsPrefix(PC,const char*[]); 200 201 extern PetscErrorCode PCComputeExplicitOperator(PC,Mat*); 202 203 /* 204 These are used to provide extra scaling of preconditioned 205 operator for time-stepping schemes like in SUNDIALS 206 */ 207 extern PetscErrorCode PCGetDiagonalScale(PC,PetscBool *); 208 extern PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec); 209 extern PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec); 210 extern PetscErrorCode PCSetDiagonalScale(PC,Vec); 211 212 /* ------------- options specific to particular preconditioners --------- */ 213 214 extern PetscErrorCode PCJacobiSetUseRowMax(PC); 215 extern PetscErrorCode PCJacobiSetUseRowSum(PC); 216 extern PetscErrorCode PCJacobiSetUseAbs(PC); 217 extern PetscErrorCode PCSORSetSymmetric(PC,MatSORType); 218 extern PetscErrorCode PCSORSetOmega(PC,PetscReal); 219 extern PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt); 220 221 extern PetscErrorCode PCEisenstatSetOmega(PC,PetscReal); 222 extern PetscErrorCode PCEisenstatNoDiagonalScaling(PC); 223 224 #define USE_PRECONDITIONER_MATRIX 0 225 #define USE_TRUE_MATRIX 1 226 extern PetscErrorCode PCBJacobiSetUseTrueLocal(PC); 227 extern PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]); 228 extern PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]); 229 230 extern PetscErrorCode PCKSPSetUseTrue(PC); 231 232 extern PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec)); 233 extern PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)); 234 extern PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec)); 235 extern PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC)); 236 extern PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*)); 237 extern PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer)); 238 extern PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC)); 239 extern PetscErrorCode PCShellGetContext(PC,void**); 240 extern PetscErrorCode PCShellSetContext(PC,void*); 241 extern PetscErrorCode PCShellSetName(PC,const char[]); 242 extern PetscErrorCode PCShellGetName(PC,char*[]); 243 244 extern PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal); 245 246 extern PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType); 247 extern PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal); 248 249 extern PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage); 250 extern PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*); 251 extern PetscErrorCode PCFactorSetUpMatSolverPackage(PC); 252 253 extern PetscErrorCode PCFactorSetFill(PC,PetscReal); 254 extern PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal); 255 extern PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal); 256 257 extern PetscErrorCode PCFactorSetMatOrderingType(PC,const MatOrderingType); 258 extern PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool ); 259 extern PetscErrorCode PCFactorSetReuseFill(PC,PetscBool ); 260 extern PetscErrorCode PCFactorSetUseInPlace(PC); 261 extern PetscErrorCode PCFactorSetAllowDiagonalFill(PC); 262 extern PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool ); 263 264 extern PetscErrorCode PCFactorSetLevels(PC,PetscInt); 265 extern PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt); 266 267 extern PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]); 268 extern PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]); 269 extern PetscErrorCode PCASMSetOverlap(PC,PetscInt); 270 extern PetscErrorCode PCASMSetSortIndices(PC,PetscBool ); 271 272 /*E 273 PCASMType - Type of additive Schwarz method to use 274 275 $ PC_ASM_BASIC - symmetric version where residuals from the ghost points are used 276 $ and computed values in ghost regions are added together. Classical 277 $ standard additive Schwarz 278 $ PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost 279 $ region are discarded. Default 280 $ PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost 281 $ region are added back in 282 $ PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded 283 $ not very good. 284 285 Level: beginner 286 287 .seealso: PCASMSetType() 288 E*/ 289 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 290 extern const char *PCASMTypes[]; 291 292 extern PetscErrorCode PCASMSetType(PC,PCASMType); 293 extern PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]); 294 extern PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]); 295 extern PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 296 extern PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]); 297 extern PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]); 298 299 /*E 300 PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per domain) 301 302 $ PC_GASM_BASIC - symmetric version where residuals from the ghost points are used 303 $ and computed values in ghost regions are added together. Classical 304 $ standard additive Schwarz 305 $ PC_GASM_RESTRICT - residuals from ghost points are used but computed values in ghost 306 $ region are discarded. Default 307 $ PC_GASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost 308 $ region are added back in 309 $ PC_GASM_NONE - ghost point residuals are not used, computed ghost values are discarded 310 $ not very good. 311 312 Level: beginner 313 314 .seealso: PCGASMSetType() 315 E*/ 316 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType; 317 extern const char *PCGASMTypes[]; 318 319 extern PetscErrorCode PCGASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]); 320 extern PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt); 321 extern PetscErrorCode PCGASMSetOverlap(PC,PetscInt); 322 extern PetscErrorCode PCGASMSetSortIndices(PC,PetscBool ); 323 324 extern PetscErrorCode PCGASMSetType(PC,PCGASMType); 325 extern PetscErrorCode PCGASMCreateSubdomains(Mat,PetscInt,IS*[]); 326 extern PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]); 327 extern PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 328 extern PetscErrorCode PCGASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]); 329 extern PetscErrorCode PCGASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]); 330 331 /*E 332 PCCompositeType - Determines how two or more preconditioner are composed 333 334 $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 335 $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 336 $ computed after the previous preconditioner application 337 $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 338 $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions) 339 $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 340 $ where first preconditioner is built from alpha I + S and second from 341 $ alpha I + R 342 343 Level: beginner 344 345 .seealso: PCCompositeSetType() 346 E*/ 347 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType; 348 extern const char *PCCompositeTypes[]; 349 350 extern PetscErrorCode PCCompositeSetUseTrue(PC); 351 extern PetscErrorCode PCCompositeSetType(PC,PCCompositeType); 352 extern PetscErrorCode PCCompositeAddPC(PC,PCType); 353 extern PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *); 354 extern PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar); 355 356 extern PetscErrorCode PCRedundantSetNumber(PC,PetscInt); 357 extern PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter); 358 extern PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*); 359 360 extern PetscErrorCode PCSPAISetEpsilon(PC,double); 361 extern PetscErrorCode PCSPAISetNBSteps(PC,PetscInt); 362 extern PetscErrorCode PCSPAISetMax(PC,PetscInt); 363 extern PetscErrorCode PCSPAISetMaxNew(PC,PetscInt); 364 extern PetscErrorCode PCSPAISetBlockSize(PC,PetscInt); 365 extern PetscErrorCode PCSPAISetCacheSize(PC,PetscInt); 366 extern PetscErrorCode PCSPAISetVerbose(PC,PetscInt); 367 extern PetscErrorCode PCSPAISetSp(PC,PetscInt); 368 369 extern PetscErrorCode PCHYPRESetType(PC,const char[]); 370 extern PetscErrorCode PCHYPREGetType(PC,const char*[]); 371 extern PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]); 372 extern PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]); 373 374 extern PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*); 375 extern PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType); 376 extern PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt); 377 extern PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS); 378 extern PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*); 379 380 /*E 381 PCFieldSplitSchurPreType - Determines how to precondition Schur complement 382 383 Level: intermediate 384 385 .seealso: PCFieldSplitSchurPrecondition() 386 E*/ 387 typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_DIAG,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType; 388 extern const char *const PCFieldSplitSchurPreTypes[]; 389 390 extern PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat); 391 extern PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*); 392 393 extern PetscErrorCode PCGalerkinSetRestriction(PC,Mat); 394 extern PetscErrorCode PCGalerkinSetInterpolation(PC,Mat); 395 396 extern PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscReal*); 397 extern PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *); 398 399 extern PetscErrorCode PCPythonSetType(PC,const char[]); 400 401 extern PetscErrorCode PCSetDM(PC,DM); 402 extern PetscErrorCode PCGetDM(PC,DM*); 403 404 extern PetscErrorCode PCSetApplicationContext(PC,void*); 405 extern PetscErrorCode PCGetApplicationContext(PC,void*); 406 407 extern PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal); 408 extern PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt); 409 extern PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool); 410 411 extern PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal); 412 extern PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool); 413 extern PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt); 414 extern PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt); 415 /*E 416 PCPARMSGlobalType - Determines the global preconditioner method in PARMS 417 418 Level: intermediate 419 420 .seealso: PCPARMSSetGlobal() 421 E*/ 422 typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType; 423 extern const char *PCPARMSGlobalTypes[]; 424 /*E 425 PCPARMSLocalType - Determines the local preconditioner method in PARMS 426 427 Level: intermediate 428 429 .seealso: PCPARMSSetLocal() 430 E*/ 431 typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType; 432 extern const char *PCPARMSLocalTypes[]; 433 434 extern PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type); 435 extern PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type); 436 extern PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits); 437 extern PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart); 438 extern PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym); 439 extern PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2); 440 441 PETSC_EXTERN_CXX_END 442 443 #endif /* __PETSCPC_H */ 444