1d03aef70SBarry Smith /* 237f753daSBarry Smith Preconditioner module. 3d03aef70SBarry Smith */ 40a835dfdSSatish Balay #if !defined(__PETSCPC_H) 50a835dfdSSatish Balay #define __PETSCPC_H 61e25c274SJed Brown #include <petscmat.h> 71e25c274SJed Brown #include <petscdmtypes.h> 8d03aef70SBarry Smith 9607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PCInitializePackage(void); 101dbb0a54SBarry Smith 11eec0b4cfSBarry Smith /* 12eec0b4cfSBarry Smith PCList contains the list of preconditioners currently registered 13bdf89e91SBarry Smith These are added with PCRegister() 14eec0b4cfSBarry Smith */ 15140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList PCList; 1682bf6240SBarry Smith 173d957683SBarry Smith /*S 188f6c3df8SBarry Smith PC - Abstract PETSc object that manages all preconditioners including direct solvers such as PCLU 193d957683SBarry Smith 203d957683SBarry Smith Level: beginner 213d957683SBarry Smith 223d957683SBarry Smith Concepts: preconditioners 233d957683SBarry Smith 241a480d89SAdministrator .seealso: PCCreate(), PCSetType(), PCType (for list of available types) 253d957683SBarry Smith S*/ 263d957683SBarry Smith typedef struct _p_PC* PC; 273d957683SBarry Smith 2876bdecfbSBarry Smith /*J 298f6c3df8SBarry Smith PCType - String with the name of a PETSc preconditioner method. 303d957683SBarry Smith 313d957683SBarry Smith Level: beginner 323d957683SBarry Smith 331a480d89SAdministrator Notes: Click on the links below to see details on a particular solver 341a480d89SAdministrator 358f6c3df8SBarry Smith PCRegister() is used to register preconditioners that are then accessible via PCSetType() 368f6c3df8SBarry Smith 378f6c3df8SBarry Smith .seealso: PCSetType(), PC, PCCreate(), PCRegister(), PCSetFromOptions() 3876bdecfbSBarry Smith J*/ 3919fd82e9SBarry Smith typedef const char* PCType; 4082bf6240SBarry Smith #define PCNONE "none" 4182bf6240SBarry Smith #define PCJACOBI "jacobi" 4282bf6240SBarry Smith #define PCSOR "sor" 4382bf6240SBarry Smith #define PCLU "lu" 4482bf6240SBarry Smith #define PCSHELL "shell" 4582bf6240SBarry Smith #define PCBJACOBI "bjacobi" 4682bf6240SBarry Smith #define PCMG "mg" 4782bf6240SBarry Smith #define PCEISENSTAT "eisenstat" 4882bf6240SBarry Smith #define PCILU "ilu" 4982bf6240SBarry Smith #define PCICC "icc" 5082bf6240SBarry Smith #define PCASM "asm" 51ab3e923fSDmitry Karpeev #define PCGASM "gasm" 5294b7f48cSBarry Smith #define PCKSP "ksp" 5382bf6240SBarry Smith #define PCCOMPOSITE "composite" 54421c37bdSBarry Smith #define PCREDUNDANT "redundant" 5527b520f0SBarry Smith #define PCSPAI "spai" 56186905e3SBarry Smith #define PCNN "nn" 574bbc92c1SBarry Smith #define PCCHOLESKY "cholesky" 583a7fca6bSBarry Smith #define PCPBJACOBI "pbjacobi" 597f5ff6fdSBarry Smith #define PCMAT "mat" 60c4888f26SBarry Smith #define PCHYPRE "hypre" 6137f80224SJose E. Roman #define PCPARMS "parms" 620971522cSBarry Smith #define PCFIELDSPLIT "fieldsplit" 63be16f70fSSatish Balay #define PCTFS "tfs" 645582bec1SHong Zhang #define PCML "ml" 652a6744ebSBarry Smith #define PCGALERKIN "galerkin" 667233f9f0SBarry Smith #define PCEXOTIC "exotic" 67bad7cb1dSBarry Smith #define PCHMPI "hmpi" 68cf037197Ssdaitch #define PCSUPPORTGRAPH "supportgraph" 69f4b8409dSBarry Smith #define PCASA "asa" 7024c02a0fSBarry Smith #define PCCP "cp" 71628b8437SMatthew Knepley #define PCBFBT "bfbt" 72519d70e2SJed Brown #define PCLSC "lsc" 731d6018f0SLisandro Dalcin #define PCPYTHON "python" 74f91d8e95SBarry Smith #define PCPFMG "pfmg" 75d851a50bSGlenn Hammond #define PCSYSPFMG "syspfmg" 76df826632SBarry Smith #define PCREDISTRIBUTE "redistribute" 77ac793be5SBarry Smith #define PCSVD "svd" 78ac793be5SBarry Smith #define PCGAMG "gamg" 79ac793be5SBarry Smith #define PCSACUSP "sacusp" /* these four run on NVIDIA GPUs using CUSP */ 80bfc29b71SVictor Minden #define PCSACUSPPOLY "sacusppoly" 818154be41SBarry Smith #define PCBICGSTABCUSP "bicgstabcusp" 8204b59e88SVictor Minden #define PCAINVCUSP "ainvcusp" 830c7d97c5SJed Brown #define PCBDDC "bddc" 84123ea438SMatthew Knepley 85123ea438SMatthew Knepley /* Logging support */ 86014dd563SJed Brown PETSC_EXTERN PetscClassId PC_CLASSID; 87123ea438SMatthew Knepley 883d957683SBarry Smith /*E 893d957683SBarry Smith PCSide - If the preconditioner is to be applied to the left, right 903d957683SBarry Smith or symmetrically around the operator. 91d03aef70SBarry Smith 923d957683SBarry Smith Level: beginner 9321e95762SBarry Smith 943d957683SBarry Smith .seealso: 953d957683SBarry Smith E*/ 969e568555SJed Brown typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide; 979e568555SJed Brown #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 980910f4e4SJed Brown PETSC_EXTERN const char *const *const PCSides; 9972b7852fSLois Curfman McInnes 100014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*); 10119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetUp(PC); 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC); 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec); 106014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec); 107014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec); 108014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec); 109014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *); 110014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec); 1114d0a8057SBarry Smith 11255849f57SBarry Smith #define PC_FILE_CLASSID 1211222 11355849f57SBarry Smith 1144d0a8057SBarry Smith /*E 1154d0a8057SBarry Smith PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates 1164d0a8057SBarry Smith 1174d0a8057SBarry Smith Level: advanced 1184d0a8057SBarry Smith 1194d0a8057SBarry Smith Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h 1204d0a8057SBarry Smith 1214d0a8057SBarry Smith .seealso: PCApplyRichardson() 1224d0a8057SBarry Smith E*/ 1234d0a8057SBarry Smith typedef enum { 1244d0a8057SBarry Smith PCRICHARDSON_CONVERGED_RTOL = 2, 1254d0a8057SBarry Smith PCRICHARDSON_CONVERGED_ATOL = 3, 1264d0a8057SBarry Smith PCRICHARDSON_CONVERGED_ITS = 4, 1274d0a8057SBarry Smith PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason; 1284d0a8057SBarry Smith 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*); 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *); 131014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool ); 13249517cdeSBarry Smith PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool); 13349517cdeSBarry Smith PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*); 13484cb2905SBarry Smith 135607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PCRegisterAll(void); 136014dd563SJed Brown PETSC_EXTERN PetscBool PCRegisterAllCalled; 13784cb2905SBarry Smith 138bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode PCRegister(const char[],PetscErrorCode(*)(PC)); 13930de9b25SBarry Smith 140014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCReset(PC); 141014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDestroy(PC*); 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC); 14319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*); 14414c91fddSBarry Smith 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*); 146014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*); 147014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*); 1485b116368SBarry Smith 149014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure); 150014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*); 151014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *); 1524b0e389bSBarry Smith 153014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer); 15455849f57SBarry Smith PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer); 1557bc3d0afSSatish Balay 156014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]); 157014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]); 158014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]); 1598ed539a5SBarry Smith 160014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*); 16171601f6fSBarry Smith 162d6913704SBarry Smith /* 163d6913704SBarry Smith These are used to provide extra scaling of preconditioned 1640f3b3ca1SHong Zhang operator for time-stepping schemes like in SUNDIALS 165d6913704SBarry Smith */ 166014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *); 167014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec); 168014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec); 169014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec); 170d6913704SBarry Smith 17184cb2905SBarry Smith /* ------------- options specific to particular preconditioners --------- */ 172329f5518SBarry Smith 173014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC); 174014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC); 175014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC); 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType); 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal); 178014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt); 179d03aef70SBarry Smith 180014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal); 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC); 182421c37bdSBarry Smith 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]); 184014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]); 1851eb62cbbSBarry Smith 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec)); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)); 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec)); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC)); 190014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*)); 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer)); 192014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC)); 193014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**); 194014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*); 195014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]); 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]); 197aabeff55SBarry Smith 198014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal); 199d90ac83dSHong Zhang 200014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType); 201014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal); 202d90ac83dSHong Zhang 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC); 2062401956bSBarry Smith 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal); 208014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal); 209014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal); 210421c37bdSBarry Smith 21119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType); 212014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool ); 213014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool ); 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC); 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC); 216014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool ); 217f5a88c2aSLois Curfman McInnes 2182591b318SToby Isaac PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC,PetscInt*); 219014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt); 220014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt); 221b35a507dSBarry Smith 222014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]); 223014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]); 224014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt); 225d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool); 226d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*); 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool ); 228f746d493SDmitry Karpeev 2293d957683SBarry Smith /*E 2303d957683SBarry Smith PCASMType - Type of additive Schwarz method to use 2313d957683SBarry Smith 232feb221f8SDmitry Karpeev $ PC_ASM_BASIC - Symmetric version where residuals from the ghost points are used 233feb221f8SDmitry Karpeev $ and computed values in ghost regions are added together. 234feb221f8SDmitry Karpeev $ Classical standard additive Schwarz. 235feb221f8SDmitry Karpeev $ PC_ASM_RESTRICT - Residuals from ghost points are used but computed values in ghost 236feb221f8SDmitry Karpeev $ region are discarded. 237feb221f8SDmitry Karpeev $ Default. 238feb221f8SDmitry Karpeev $ PC_ASM_INTERPOLATE - Residuals from ghost points are not used, computed values in ghost 239feb221f8SDmitry Karpeev $ region are added back in. 240feb221f8SDmitry Karpeev $ PC_ASM_NONE - Residuals from ghost points are not used, computed ghost values are 241feb221f8SDmitry Karpeev $ discarded. 242feb221f8SDmitry Karpeev $ Not very good. 2433d957683SBarry Smith 2443d957683SBarry Smith Level: beginner 2453d957683SBarry Smith 2463d957683SBarry Smith .seealso: PCASMSetType() 2473d957683SBarry Smith E*/ 248d252947aSBarry Smith typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 2496a6fc655SJed Brown PETSC_EXTERN const char *const PCASMTypes[]; 2509dcbbd2bSBarry Smith 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType); 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]); 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]); 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]); 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]); 257981c4779SBarry Smith 2583d957683SBarry Smith /*E 2596a4f0f73SDmitry Karpeev PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain). 260f746d493SDmitry Karpeev 2616a4f0f73SDmitry Karpeev Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 2626a4f0f73SDmitry Karpeev domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 2636a4f0f73SDmitry Karpeev a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 2646a4f0f73SDmitry Karpeev (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 2656a4f0f73SDmitry Karpeev 2666a4f0f73SDmitry Karpeev $ PC_GASM_BASIC - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 2676a4f0f73SDmitry Karpeev $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 2686a4f0f73SDmitry Karpeev $ from neighboring subdomains. 269feb221f8SDmitry Karpeev $ Classical standard additive Schwarz. 2706a4f0f73SDmitry Karpeev $ PC_GASM_RESTRICT - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 2716a4f0f73SDmitry Karpeev $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 2726a4f0f73SDmitry Karpeev $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 2736a4f0f73SDmitry Karpeev $ assumption). 274feb221f8SDmitry Karpeev $ Default. 2756a4f0f73SDmitry Karpeev $ PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 2766a4f0f73SDmitry Karpeev $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 2776a4f0f73SDmitry Karpeev $ from neighboring subdomains. 2786a4f0f73SDmitry Karpeev $ 2796a4f0f73SDmitry Karpeev $ PC_GASM_NONE - Residuals and corrections are zeroed out outside the local subdomains. 280feb221f8SDmitry Karpeev $ Not very good. 281f746d493SDmitry Karpeev 282f746d493SDmitry Karpeev Level: beginner 283f746d493SDmitry Karpeev 284f746d493SDmitry Karpeev .seealso: PCGASMSetType() 285f746d493SDmitry Karpeev E*/ 286f746d493SDmitry Karpeev typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType; 2876a6fc655SJed Brown PETSC_EXTERN const char *const PCGASMTypes[]; 288f746d493SDmitry Karpeev 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool); 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt); 292d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCGASMSetDMSubdomains(PC,PetscBool); 293d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCGASMGetDMSubdomains(PC,PetscBool*); 294014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool ); 295f746d493SDmitry Karpeev 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType); 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]); 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]); 299014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 300014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]); 301014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]); 302f746d493SDmitry Karpeev 303f746d493SDmitry Karpeev /*E 3043d957683SBarry Smith PCCompositeType - Determines how two or more preconditioner are composed 3053d957683SBarry Smith 3063d957683SBarry Smith $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 3073d957683SBarry Smith $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 3083d957683SBarry Smith $ computed after the previous preconditioner application 309ccb205f8SBarry Smith $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 310ccb205f8SBarry Smith $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions) 3113d957683SBarry Smith $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 3123d957683SBarry Smith $ where first preconditioner is built from alpha I + S and second from 3133d957683SBarry Smith $ alpha I + R 3143d957683SBarry Smith 3153d957683SBarry Smith Level: beginner 3163d957683SBarry Smith 3173d957683SBarry Smith .seealso: PCCompositeSetType() 3183d957683SBarry Smith E*/ 3193b224e63SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType; 3206a6fc655SJed Brown PETSC_EXTERN const char *const PCCompositeTypes[]; 3219dcbbd2bSBarry Smith 322014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType); 32319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType); 324014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *); 325014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar); 326981c4779SBarry Smith 327014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt); 328014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter); 329014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*); 330da3a660dSBarry Smith 331014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double); 332014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt); 333014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt); 334014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt); 335014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt); 336014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt); 337014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt); 338014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt); 3393304466cSBarry Smith 340014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]); 341014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]); 342014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]); 343014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]); 3443304466cSBarry Smith 345014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*); 346014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*); 347014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType); 348014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt); 349014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS); 350014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*); 3514ab8060aSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool); 3524ab8060aSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*); 353084e4875SJed Brown 354084e4875SJed Brown /*E 355084e4875SJed Brown PCFieldSplitSchurPreType - Determines how to precondition Schur complement 356084e4875SJed Brown 357084e4875SJed Brown Level: intermediate 358084e4875SJed Brown 359084e4875SJed Brown .seealso: PCFieldSplitSchurPrecondition() 360084e4875SJed Brown E*/ 361e87fab1bSBarry Smith typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_A11,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType; 362014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[]; 363084e4875SJed Brown 364ab1df9f5SJed Brown /*E 365c9c6ffaaSJed Brown PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 366ab1df9f5SJed Brown 367ab1df9f5SJed Brown Level: intermediate 368ab1df9f5SJed Brown 369c9c6ffaaSJed Brown .seealso: PCFieldSplitSetSchurFactType() 370ab1df9f5SJed Brown E*/ 371ab1df9f5SJed Brown typedef enum { 372c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_DIAG, 373c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_LOWER, 374c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_UPPER, 375c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_FULL 376c9c6ffaaSJed Brown } PCFieldSplitSchurFactType; 377014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[]; 378ab1df9f5SJed Brown 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*); 3823d30b1ffSBarry Smith 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat); 3852a6744ebSBarry Smith 386014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*); 387014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *); 3882102ba4dSSatish Balay 389014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]); 39067fac13cSBarry Smith 391014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM); 392014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*); 3936c699258SBarry Smith 394014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*); 395014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*); 3961b2093e4SBarry Smith 397014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal); 398014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt); 399014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool); 4006c699258SBarry Smith 401014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal); 402014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool); 403014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt); 404014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt); 40537f80224SJose E. Roman /*E 40637f80224SJose E. Roman PCPARMSGlobalType - Determines the global preconditioner method in PARMS 40737f80224SJose E. Roman 40837f80224SJose E. Roman Level: intermediate 40937f80224SJose E. Roman 41037f80224SJose E. Roman .seealso: PCPARMSSetGlobal() 41137f80224SJose E. Roman E*/ 41237f80224SJose E. Roman typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType; 4136a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSGlobalTypes[]; 41437f80224SJose E. Roman /*E 41537f80224SJose E. Roman PCPARMSLocalType - Determines the local preconditioner method in PARMS 41637f80224SJose E. Roman 41737f80224SJose E. Roman Level: intermediate 41837f80224SJose E. Roman 41937f80224SJose E. Roman .seealso: PCPARMSSetLocal() 42037f80224SJose E. Roman E*/ 42137f80224SJose E. Roman typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType; 4226a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSLocalTypes[]; 42337f80224SJose E. Roman 424014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type); 425014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type); 426014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits); 427014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart); 428014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym); 429014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2); 43037f80224SJose E. Roman 431*5154939eSJed Brown /*E 432*5154939eSJed Brown PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method 433*5154939eSJed Brown 434*5154939eSJed Brown Level: intermediate 435*5154939eSJed Brown 436*5154939eSJed Brown .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseProl() 437*5154939eSJed Brown E*/ 438a782d3a5SMark Adams typedef const char *PCGAMGType; 439a782d3a5SMark Adams #define PCGAMGAGG "agg" 440a782d3a5SMark Adams #define PCGAMGGEO "geo" 441014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt); 442014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool); 443014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool); 444014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt); 445014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal); 446014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt); 447014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt); 44819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType ); 449014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n); 450014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n); 451014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool); 452fc897844SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool); 4533e3471ccSMark Adams PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void); 4543e3471ccSMark Adams PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void); 45556de70efSSatish Balay 4560c7d97c5SJed Brown #if defined(PETSC_HAVE_PCBDDC) 45753cdbc3dSStefano Zampini /* Enum defining how to treat the coarse problem */ 4580c7d97c5SJed Brown typedef enum {SEQUENTIAL_BDDC,REPLICATED_BDDC,PARALLEL_BDDC,MULTILEVEL_BDDC} CoarseProblemType; 4594fad6a16SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt); 4604fad6a16SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetMaxLevels(PC,PetscInt); 4610bdf917eSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace); 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS); 463014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*); 464014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS); 465014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*); 466014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseProblemType(PC,CoarseProblemType); 467014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]); 4681a83f524SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode); 4693425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*); 4703425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec); 4713425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec); 4720c7d97c5SJed Brown #endif 4730c7d97c5SJed Brown 474b965d45eSStefano Zampini PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool); 475014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar); 47646caae47SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec); 477d3b1e0e7SMatthew Knepley 478b4876bcbSBarry Smith /*E 479b4876bcbSBarry Smith PCMGType - Determines the type of multigrid method that is run. 480b4876bcbSBarry Smith 481b4876bcbSBarry Smith Level: beginner 482b4876bcbSBarry Smith 483b4876bcbSBarry Smith Values: 484b4876bcbSBarry Smith + PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles() 485b4876bcbSBarry Smith . PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are 486b4876bcbSBarry Smith smoothed before updating the residual. This only uses the 487b4876bcbSBarry Smith down smoother, in the preconditioner the upper smoother is ignored 488b4876bcbSBarry Smith . PC_MG_FULL - same as multiplicative except one also performs grid sequencing, 489b4876bcbSBarry Smith that is starts on the coarsest grid, performs a cycle, interpolates 490b4876bcbSBarry Smith to the next, performs a cycle etc. This is much like the F-cycle presented in "Multigrid" by Trottenberg, Oosterlee, Schuller page 49, but that 491b4876bcbSBarry Smith algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 492b4876bcbSBarry Smith calls the V-cycle only on the coarser level and has a post-smoother instead. 493b4876bcbSBarry Smith - PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level 494b4876bcbSBarry Smith from a finer 495b4876bcbSBarry Smith 496b4876bcbSBarry Smith .seealso: PCMGSetType() 497b4876bcbSBarry Smith 498b4876bcbSBarry Smith E*/ 499b4876bcbSBarry Smith typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType; 500b4876bcbSBarry Smith PETSC_EXTERN const char *const PCMGTypes[]; 501b4876bcbSBarry Smith #define PC_MG_CASCADE PC_MG_KASKADE; 502b4876bcbSBarry Smith 503b4876bcbSBarry Smith /*E 504b4876bcbSBarry Smith PCMGCycleType - Use V-cycle or W-cycle 505b4876bcbSBarry Smith 506b4876bcbSBarry Smith Level: beginner 507b4876bcbSBarry Smith 508b4876bcbSBarry Smith Values: 509b4876bcbSBarry Smith + PC_MG_V_CYCLE 510b4876bcbSBarry Smith - PC_MG_W_CYCLE 511b4876bcbSBarry Smith 512b4876bcbSBarry Smith .seealso: PCMGSetCycleType() 513b4876bcbSBarry Smith 514b4876bcbSBarry Smith E*/ 515b4876bcbSBarry Smith typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType; 516b4876bcbSBarry Smith PETSC_EXTERN const char *const PCMGCycleTypes[]; 517b4876bcbSBarry Smith 518b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType); 519b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*); 520b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*); 521b4876bcbSBarry Smith 522b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt); 523b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt); 524b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType); 525b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType); 526b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt); 527b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt); 528b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool); 529b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*); 530b4876bcbSBarry Smith 531b4876bcbSBarry Smith 532b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec); 533b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec); 534b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec); 535b4876bcbSBarry Smith 536b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat); 537b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*); 538b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat); 539b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*); 540b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec); 541b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*); 542b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat); 543b4876bcbSBarry Smith 544b4876bcbSBarry Smith /*E 545b4876bcbSBarry Smith PCExoticType - Face based or wirebasket based coarse grid space 546b4876bcbSBarry Smith 547b4876bcbSBarry Smith Level: beginner 548b4876bcbSBarry Smith 549b4876bcbSBarry Smith .seealso: PCExoticSetType(), PCEXOTIC 550b4876bcbSBarry Smith E*/ 551b4876bcbSBarry Smith typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType; 552b4876bcbSBarry Smith PETSC_EXTERN const char *const PCExoticTypes[]; 553b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType); 554b4876bcbSBarry Smith 555123ea438SMatthew Knepley #endif /* __PETSCPC_H */ 556