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