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