1 /* 2 Preconditioner module. 3 */ 4 #if !defined(__PETSCPC_H) 5 #define __PETSCPC_H 6 #include "petscmat.h" 7 PETSC_EXTERN_CXX_BEGIN 8 9 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT 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 PCKSP "ksp" 52 #define PCCOMPOSITE "composite" 53 #define PCREDUNDANT "redundant" 54 #define PCSPAI "spai" 55 #define PCNN "nn" 56 #define PCCHOLESKY "cholesky" 57 #define PCSAMG "samg" 58 #define PCPBJACOBI "pbjacobi" 59 #define PCMAT "mat" 60 #define PCHYPRE "hypre" 61 #define PCFIELDSPLIT "fieldsplit" 62 #define PCTFS "tfs" 63 #define PCML "ml" 64 #define PCPROMETHEUS "prometheus" 65 #define PCGALERKIN "galerkin" 66 #define PCEXOTIC "exotic" 67 #define PCOPENMP "openmp" 68 #define PCSUPPORTGRAPH "supportgraph" 69 #define PCASA "asa" 70 #define PCCP "cp" 71 #define PCBFBT "bfbt" 72 #define PCPYTHON "python" 73 74 /* Logging support */ 75 extern PetscCookie PETSCKSP_DLLEXPORT PC_COOKIE; 76 77 /*E 78 PCSide - If the preconditioner is to be applied to the left, right 79 or symmetrically around the operator. 80 81 Level: beginner 82 83 .seealso: 84 E*/ 85 typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide; 86 extern const char *PCSides[]; 87 88 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCreate(MPI_Comm,PC*); 89 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetType(PC,const PCType); 90 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetUp(PC); 91 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetUpOnBlocks(PC); 92 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApply(PC,Vec,Vec); 93 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricLeft(PC,Vec,Vec); 94 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricRight(PC,Vec,Vec); 95 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec); 96 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyTranspose(PC,Vec,Vec); 97 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyTransposeExists(PC,PetscTruth*); 98 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec); 99 100 /*E 101 PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates 102 103 Level: advanced 104 105 Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h 106 107 .seealso: PCApplyRichardson() 108 E*/ 109 typedef enum { 110 PCRICHARDSON_CONVERGED_RTOL = 2, 111 PCRICHARDSON_CONVERGED_ATOL = 3, 112 PCRICHARDSON_CONVERGED_ITS = 4, 113 PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason; 114 115 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscInt*,PCRichardsonConvergedReason*); 116 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardsonExists(PC,PetscTruth*); 117 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetInitialGuessNonzero(PC,PetscTruth); 118 119 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterDestroy(void); 120 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterAll(const char[]); 121 extern PetscTruth PCRegisterAllCalled; 122 123 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC)); 124 125 /*MC 126 PCRegisterDynamic - Adds a method to the preconditioner package. 127 128 Synopsis: 129 PetscErrorCode PCRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(PC)) 130 131 Not collective 132 133 Input Parameters: 134 + name_solver - name of a new user-defined solver 135 . path - path (either absolute or relative) the library containing this solver 136 . name_create - name of routine to create method context 137 - routine_create - routine to create method context 138 139 Notes: 140 PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners. 141 142 If dynamic libraries are used, then the fourth input argument (routine_create) 143 is ignored. 144 145 Sample usage: 146 .vb 147 PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib", 148 "MySolverCreate",MySolverCreate); 149 .ve 150 151 Then, your solver can be chosen with the procedural interface via 152 $ PCSetType(pc,"my_solver") 153 or at runtime via the option 154 $ -pc_type my_solver 155 156 Level: advanced 157 158 Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, or ${any environmental variable} 159 occuring in pathname will be replaced with appropriate values. 160 If your function is not being put into a shared library then use PCRegister() instead 161 162 .keywords: PC, register 163 164 .seealso: PCRegisterAll(), PCRegisterDestroy() 165 M*/ 166 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 167 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0) 168 #else 169 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d) 170 #endif 171 172 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDestroy(PC); 173 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetFromOptions(PC); 174 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetType(PC,const PCType*); 175 176 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatrix(PC,Mat*); 177 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*); 178 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*); 179 180 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetOperators(PC,Mat,Mat,MatStructure); 181 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperators(PC,Mat*,Mat*,MatStructure*); 182 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperatorsSet(PC,PetscTruth*,PetscTruth*); 183 184 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCView(PC,PetscViewer); 185 186 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetOptionsPrefix(PC,const char[]); 187 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCAppendOptionsPrefix(PC,const char[]); 188 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOptionsPrefix(PC,const char*[]); 189 190 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCComputeExplicitOperator(PC,Mat*); 191 192 /* 193 These are used to provide extra scaling of preconditioned 194 operator for time-stepping schemes like in SUNDIALS 195 */ 196 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScale(PC,PetscTruth*); 197 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleLeft(PC,Vec,Vec); 198 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleRight(PC,Vec,Vec); 199 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleSet(PC,Vec); 200 201 /* ------------- options specific to particular preconditioners --------- */ 202 203 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowMax(PC); 204 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowSum(PC); 205 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseAbs(PC); 206 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric(PC,MatSORType); 207 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega(PC,PetscReal); 208 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations(PC,PetscInt,PetscInt); 209 210 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatSetOmega(PC,PetscReal); 211 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatNoDiagonalScaling(PC); 212 213 #define USE_PRECONDITIONER_MATRIX 0 214 #define USE_TRUE_MATRIX 1 215 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetUseTrueLocal(PC); 216 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]); 217 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]); 218 219 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPSetUseTrue(PC); 220 221 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApply(PC,PetscErrorCode (*)(void*,Vec,Vec)); 222 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyBA(PC,PetscErrorCode (*)(void*,PCSide,Vec,Vec,Vec)); 223 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyTranspose(PC,PetscErrorCode (*)(void*,Vec,Vec)); 224 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetSetUp(PC,PetscErrorCode (*)(void*)); 225 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyRichardson(PC,PetscErrorCode (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscInt*,PCRichardsonConvergedReason*)); 226 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetView(PC,PetscErrorCode (*)(void*,PetscViewer)); 227 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetDestroy(PC,PetscErrorCode (*)(void*)); 228 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetContext(PC,void**); 229 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetContext(PC,void*); 230 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetName(PC,const char[]); 231 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetName(PC,char*[]); 232 233 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot(PC,PetscReal); 234 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero(PC,PetscReal); 235 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd(PC,PetscTruth); 236 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC,const MatSolverPackage); 237 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage(PC,const MatSolverPackage*); 238 239 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC,PetscReal); 240 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC,PetscReal); 241 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC,PetscReal); 242 243 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC,const MatOrderingType); 244 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC,PetscTruth); 245 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC,PetscTruth); 246 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC); 247 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC); 248 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC,PetscTruth); 249 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks(PC,PetscReal); 250 251 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC,PetscInt); 252 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC,PetscReal,PetscReal,PetscInt); 253 254 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetLocalSubdomains(PC,PetscInt,IS[]); 255 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetTotalSubdomains(PC,PetscInt,IS[]); 256 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetOverlap(PC,PetscInt); 257 /*E 258 PCASMType - Type of additive Schwarz method to use 259 260 $ PC_ASM_BASIC - symmetric version where residuals from the ghost points are used 261 $ and computed values in ghost regions are added together. Classical 262 $ standard additive Schwarz 263 $ PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost 264 $ region are discarded. Default 265 $ PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost 266 $ region are added back in 267 $ PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded 268 $ not very good. 269 270 Level: beginner 271 272 .seealso: PCASMSetType() 273 E*/ 274 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 275 extern const char *PCASMTypes[]; 276 277 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetType(PC,PCASMType); 278 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains(Mat,PetscInt,IS*[]); 279 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMDestroySubdomains(PetscInt,IS[]); 280 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt *,IS **); 281 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetUseInPlace(PC); 282 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubdomains(PC,PetscInt*,IS*[]); 283 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]); 284 285 /*E 286 PCCompositeType - Determines how two or more preconditioner are composed 287 288 $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 289 $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 290 $ computed after the previous preconditioner application 291 $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 292 $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions) 293 $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 294 $ where first preconditioner is built from alpha I + S and second from 295 $ alpha I + R 296 297 Level: beginner 298 299 .seealso: PCCompositeSetType() 300 E*/ 301 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType; 302 extern const char *PCCompositeTypes[]; 303 304 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC); 305 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC,PCCompositeType); 306 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC,PCType); 307 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC,PetscInt,PC *); 308 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC,PetscScalar); 309 310 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber(PC,PetscInt); 311 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter(PC,VecScatter,VecScatter); 312 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators(PC,Mat*,Mat*); 313 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC(PC,PC*); 314 315 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetEpsilon(PC,double); 316 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetNBSteps(PC,PetscInt); 317 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMax(PC,PetscInt); 318 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMaxNew(PC,PetscInt); 319 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetBlockSize(PC,PetscInt); 320 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetCacheSize(PC,PetscInt); 321 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetVerbose(PC,PetscInt); 322 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetSp(PC,PetscInt); 323 324 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC,const char[]); 325 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC,const char*[]); 326 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]); 327 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]); 328 329 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC,PetscInt,PetscInt*); 330 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC,PCCompositeType); 331 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC,PetscInt); 332 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC,IS); 333 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC,PetscTruth); 334 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*); 335 336 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetRestriction(PC,Mat); 337 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetInterpolation(PC,Mat); 338 339 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetCoordinates(PC,PetscInt,PetscReal*); 340 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSASetVectors(PC,PetscInt,PetscReal *); 341 342 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPythonSetType(PC,const char[]); 343 344 PETSC_EXTERN_CXX_END 345 #endif /* __PETSCPC_H */ 346