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 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 #define PCType char* 17 18 /*S 19 PC - Abstract PETSc object that manages all preconditioners 20 21 Level: beginner 22 23 Concepts: preconditioners 24 25 .seealso: PCCreate(), PCSetType(), PCType (for list of available types) 26 S*/ 27 typedef struct _p_PC* PC; 28 29 /*E 30 PCType - String with the name of a PETSc preconditioner method or the creation function 31 with an optional dynamic library name, for example 32 http://www.mcs.anl.gov/petsc/lib.a:mypccreate() 33 34 Level: beginner 35 36 Notes: Click on the links below to see details on a particular solver 37 38 .seealso: PCSetType(), PC, PCCreate() 39 E*/ 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 65 /* Logging support */ 66 extern PetscCookie PC_COOKIE; 67 extern PetscEvent PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft; 68 extern PetscEvent PC_ApplySymmetricRight, PC_ModifySubMatrices; 69 70 /*E 71 PCSide - If the preconditioner is to be applied to the left, right 72 or symmetrically around the operator. 73 74 Level: beginner 75 76 .seealso: 77 E*/ 78 typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide; 79 80 EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*); 81 EXTERN PetscErrorCode PCSetType(PC,const PCType); 82 EXTERN PetscErrorCode PCSetUp(PC); 83 EXTERN PetscErrorCode PCSetUpOnBlocks(PC); 84 EXTERN PetscErrorCode PCApply(PC,Vec,Vec); 85 EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec); 86 EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec); 87 EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec); 88 EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec); 89 EXTERN PetscErrorCode PCHasApplyTranspose(PC,PetscTruth*); 90 EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec); 91 EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt); 92 EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscTruth*); 93 94 EXTERN PetscErrorCode PCRegisterDestroy(void); 95 EXTERN PetscErrorCode PCRegisterAll(const char[]); 96 extern PetscTruth PCRegisterAllCalled; 97 98 EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC)); 99 100 /*MC 101 PCRegisterDynamic - Adds a method to the preconditioner package. 102 103 Synopsis: 104 PetscErrorCode PCRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(PC)) 105 106 Not collective 107 108 Input Parameters: 109 + name_solver - name of a new user-defined solver 110 . path - path (either absolute or relative) the library containing this solver 111 . name_create - name of routine to create method context 112 - routine_create - routine to create method context 113 114 Notes: 115 PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners. 116 117 If dynamic libraries are used, then the fourth input argument (routine_create) 118 is ignored. 119 120 Sample usage: 121 .vb 122 PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib", 123 "MySolverCreate",MySolverCreate); 124 .ve 125 126 Then, your solver can be chosen with the procedural interface via 127 $ PCSetType(pc,"my_solver") 128 or at runtime via the option 129 $ -pc_type my_solver 130 131 Level: advanced 132 133 Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, or ${any environmental variable} 134 occuring in pathname will be replaced with appropriate values. 135 If your function is not being put into a shared library then use PCRegister() instead 136 137 .keywords: PC, register 138 139 .seealso: PCRegisterAll(), PCRegisterDestroy() 140 M*/ 141 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 142 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0) 143 #else 144 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d) 145 #endif 146 147 EXTERN PetscErrorCode PCDestroy(PC); 148 EXTERN PetscErrorCode PCSetFromOptions(PC); 149 EXTERN PetscErrorCode PCGetType(PC,PCType*); 150 151 EXTERN PetscErrorCode PCGetFactoredMatrix(PC,Mat*); 152 EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*); 153 EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*); 154 155 EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure); 156 EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*); 157 158 EXTERN PetscErrorCode PCView(PC,PetscViewer); 159 160 EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]); 161 EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]); 162 EXTERN PetscErrorCode PCGetOptionsPrefix(PC,char*[]); 163 164 EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*); 165 166 /* 167 These are used to provide extra scaling of preconditioned 168 operator for time-stepping schemes like in PVODE 169 */ 170 EXTERN PetscErrorCode PCDiagonalScale(PC,PetscTruth*); 171 EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec); 172 EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec); 173 EXTERN PetscErrorCode PCDiagonalScaleSet(PC,Vec); 174 175 /* ------------- options specific to particular preconditioners --------- */ 176 177 EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC); 178 EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType); 179 EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal); 180 EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt); 181 182 EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal); 183 EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC); 184 185 #define USE_PRECONDITIONER_MATRIX 0 186 #define USE_TRUE_MATRIX 1 187 EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC); 188 EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]); 189 EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]); 190 191 EXTERN PetscErrorCode PCKSPSetUseTrue(PC); 192 193 EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(void*,Vec,Vec),void*); 194 EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(void*,Vec,Vec)); 195 EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(void*)); 196 EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt),void*); 197 EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(void*,PetscViewer)); 198 EXTERN PetscErrorCode PCShellSetName(PC,const char[]); 199 EXTERN PetscErrorCode PCShellGetName(PC,char*[]); 200 201 EXTERN PetscErrorCode PCLUSetMatOrdering(PC,MatOrderingType); 202 EXTERN PetscErrorCode PCLUSetReuseOrdering(PC,PetscTruth); 203 EXTERN PetscErrorCode PCLUSetReuseFill(PC,PetscTruth); 204 EXTERN PetscErrorCode PCLUSetUseInPlace(PC); 205 EXTERN PetscErrorCode PCLUSetFill(PC,PetscReal); 206 EXTERN PetscErrorCode PCLUSetDamping(PC,PetscReal); 207 EXTERN PetscErrorCode PCLUSetShift(PC,PetscTruth); 208 EXTERN PetscErrorCode PCLUSetPivoting(PC,PetscReal); 209 EXTERN PetscErrorCode PCLUSetPivotInBlocks(PC,PetscTruth); 210 EXTERN PetscErrorCode PCLUSetZeroPivot(PC,PetscReal); 211 EXTERN PetscErrorCode PCLUReorderForNonzeroDiagonal(PC,PetscReal); 212 213 EXTERN PetscErrorCode PCCholeskySetMatOrdering(PC,MatOrderingType); 214 EXTERN PetscErrorCode PCCholeskySetReuseOrdering(PC,PetscTruth); 215 EXTERN PetscErrorCode PCCholeskySetReuseFill(PC,PetscTruth); 216 EXTERN PetscErrorCode PCCholeskySetUseInPlace(PC); 217 EXTERN PetscErrorCode PCCholeskySetFill(PC,PetscReal); 218 EXTERN PetscErrorCode PCCholeskySetDamping(PC,PetscReal); 219 EXTERN PetscErrorCode PCCholeskySetShift(PC,PetscTruth); 220 EXTERN PetscErrorCode PCCholeskySetPivotInBlocks(PC,PetscTruth); 221 222 EXTERN PetscErrorCode PCILUSetMatOrdering(PC,MatOrderingType); 223 EXTERN PetscErrorCode PCILUSetUseInPlace(PC); 224 EXTERN PetscErrorCode PCILUSetFill(PC,PetscReal); 225 EXTERN PetscErrorCode PCILUSetLevels(PC,PetscInt); 226 EXTERN PetscErrorCode PCILUSetReuseOrdering(PC,PetscTruth); 227 EXTERN PetscErrorCode PCILUSetUseDropTolerance(PC,PetscReal,PetscReal,PetscInt); 228 EXTERN PetscErrorCode PCILUDTSetReuseFill(PC,PetscTruth); 229 EXTERN PetscErrorCode PCILUSetAllowDiagonalFill(PC); 230 EXTERN PetscErrorCode PCILUSetDamping(PC,PetscReal); 231 EXTERN PetscErrorCode PCFactorSetShiftNonzero(PetscReal,MatFactorInfo*); /* to be moved to other position! */ 232 EXTERN PetscErrorCode PCILUSetShift(PC,PetscTruth); 233 EXTERN PetscErrorCode PCILUSetPivotInBlocks(PC,PetscTruth); 234 EXTERN PetscErrorCode PCILUSetZeroPivot(PC,PetscReal); 235 EXTERN PetscErrorCode PCILUReorderForNonzeroDiagonal(PC,PetscReal); 236 237 EXTERN PetscErrorCode PCICCSetMatOrdering(PC,MatOrderingType); 238 EXTERN PetscErrorCode PCICCSetFill(PC,PetscReal); 239 EXTERN PetscErrorCode PCICCSetLevels(PC,PetscInt); 240 EXTERN PetscErrorCode PCICCSetDamping(PC,PetscReal); 241 EXTERN PetscErrorCode PCICCSetShift(PC,PetscTruth); 242 EXTERN PetscErrorCode PCICCSetPivotInBlocks(PC,PetscTruth); 243 EXTERN PetscErrorCode PCICCSetZeroPivot(PC,PetscReal); 244 245 EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[]); 246 EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[]); 247 EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt); 248 /*E 249 PCASMType - Type of additive Schwarz method to use 250 251 $ PC_ASM_BASIC - symmetric version where residuals from the ghost points are used 252 $ and computed values in ghost regions are added together. Classical 253 $ standard additive Schwarz 254 $ PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost 255 $ region are discarded. Default 256 $ PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost 257 $ region are added back in 258 $ PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded 259 $ not very good. 260 261 Level: beginner 262 263 .seealso: PCASMSetType() 264 E*/ 265 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 266 EXTERN PetscErrorCode PCASMSetType(PC,PCASMType); 267 EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt *,IS **); 268 EXTERN PetscErrorCode PCASMSetUseInPlace(PC); 269 EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[]); 270 EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]); 271 272 /*E 273 PCCompositeType - Determines how two or more preconditioner are composed 274 275 $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 276 $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 277 $ computed after the previous preconditioner application 278 $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 279 $ where first preconditioner is built from alpha I + S and second from 280 $ alpha I + R 281 282 Level: beginner 283 284 .seealso: PCCompositeSetType() 285 E*/ 286 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType; 287 EXTERN PetscErrorCode PCCompositeSetUseTrue(PC); 288 EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType); 289 EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType); 290 EXTERN PetscErrorCode PCCompositeGetPC(PC pc,PetscInt n,PC *); 291 EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar); 292 293 EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter); 294 EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*); 295 EXTERN PetscErrorCode PCRedundantGetPC(PC,PC*); 296 EXTERN PetscErrorCode MatGetOrderingList(PetscFList *list); 297 298 EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double); 299 EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt); 300 EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt); 301 EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt); 302 EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt); 303 EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt); 304 EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt); 305 EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt); 306 307 EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]); 308 EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]); 309 EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]); 310 311 EXTERN PetscErrorCode PCFieldSplitSetFields(PC,PetscInt,PetscInt*); 312 EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType); 313 314 PETSC_EXTERN_CXX_END 315 #endif /* __PETSCPC_H */ 316