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" 67f4b8409dSBarry Smith #define PCASA "asa" 6824c02a0fSBarry Smith #define PCCP "cp" 69628b8437SMatthew Knepley #define PCBFBT "bfbt" 70519d70e2SJed Brown #define PCLSC "lsc" 711d6018f0SLisandro Dalcin #define PCPYTHON "python" 72f91d8e95SBarry Smith #define PCPFMG "pfmg" 73d851a50bSGlenn Hammond #define PCSYSPFMG "syspfmg" 74df826632SBarry Smith #define PCREDISTRIBUTE "redistribute" 75ac793be5SBarry Smith #define PCSVD "svd" 76ac793be5SBarry Smith #define PCGAMG "gamg" 77ac793be5SBarry Smith #define PCSACUSP "sacusp" /* these four run on NVIDIA GPUs using CUSP */ 78bfc29b71SVictor Minden #define PCSACUSPPOLY "sacusppoly" 798154be41SBarry Smith #define PCBICGSTABCUSP "bicgstabcusp" 8004b59e88SVictor Minden #define PCAINVCUSP "ainvcusp" 810c7d97c5SJed Brown #define PCBDDC "bddc" 82123ea438SMatthew Knepley 83123ea438SMatthew Knepley /* Logging support */ 84014dd563SJed Brown PETSC_EXTERN PetscClassId PC_CLASSID; 85123ea438SMatthew Knepley 863d957683SBarry Smith /*E 873d957683SBarry Smith PCSide - If the preconditioner is to be applied to the left, right 883d957683SBarry Smith or symmetrically around the operator. 89d03aef70SBarry Smith 903d957683SBarry Smith Level: beginner 9121e95762SBarry Smith 923d957683SBarry Smith .seealso: 933d957683SBarry Smith E*/ 949e568555SJed Brown typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide; 959e568555SJed Brown #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 960910f4e4SJed Brown PETSC_EXTERN const char *const *const PCSides; 9772b7852fSLois Curfman McInnes 98014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*); 9919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType); 100014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetUp(PC); 101014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec); 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec); 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec); 106014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec); 107014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *); 108014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec); 1094d0a8057SBarry Smith 11055849f57SBarry Smith #define PC_FILE_CLASSID 1211222 11155849f57SBarry Smith 1124d0a8057SBarry Smith /*E 1134d0a8057SBarry Smith PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates 1144d0a8057SBarry Smith 1154d0a8057SBarry Smith Level: advanced 1164d0a8057SBarry Smith 1174d0a8057SBarry Smith Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h 1184d0a8057SBarry Smith 1194d0a8057SBarry Smith .seealso: PCApplyRichardson() 1204d0a8057SBarry Smith E*/ 1214d0a8057SBarry Smith typedef enum { 1224d0a8057SBarry Smith PCRICHARDSON_CONVERGED_RTOL = 2, 1234d0a8057SBarry Smith PCRICHARDSON_CONVERGED_ATOL = 3, 1244d0a8057SBarry Smith PCRICHARDSON_CONVERGED_ITS = 4, 1254d0a8057SBarry Smith PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason; 1264d0a8057SBarry Smith 127014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*); 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *); 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool ); 13049517cdeSBarry Smith PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool); 13149517cdeSBarry Smith PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*); 13284cb2905SBarry Smith 133607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PCRegisterAll(void); 134014dd563SJed Brown PETSC_EXTERN PetscBool PCRegisterAllCalled; 13584cb2905SBarry Smith 136bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode PCRegister(const char[],PetscErrorCode(*)(PC)); 13730de9b25SBarry Smith 138014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCReset(PC); 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDestroy(PC*); 140014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC); 14119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*); 14214c91fddSBarry Smith 143014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*); 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*); 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*); 1465b116368SBarry Smith 147014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure); 148014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*); 149014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *); 1504b0e389bSBarry Smith 151014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer); 15255849f57SBarry Smith PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer); 153ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PCViewFromOptions(PC A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 1547bc3d0afSSatish Balay 155014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]); 156014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]); 157014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]); 1588ed539a5SBarry Smith 159014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*); 16071601f6fSBarry Smith 161d6913704SBarry Smith /* 162d6913704SBarry Smith These are used to provide extra scaling of preconditioned 1630f3b3ca1SHong Zhang operator for time-stepping schemes like in SUNDIALS 164d6913704SBarry Smith */ 165014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *); 166014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec); 167014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec); 168014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec); 169d6913704SBarry Smith 17084cb2905SBarry Smith /* ------------- options specific to particular preconditioners --------- */ 171329f5518SBarry Smith 172014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC); 173014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC); 174014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC); 175014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType); 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal); 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt); 178d03aef70SBarry Smith 179014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal); 180014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC); 181421c37bdSBarry Smith 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]); 1841eb62cbbSBarry Smith 185014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec)); 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec)); 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC)); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*)); 190014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer)); 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC)); 192014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**); 193014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*); 194014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]); 195014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]); 196aabeff55SBarry Smith 197014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal); 198d90ac83dSHong Zhang 199014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType); 200014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal); 201d90ac83dSHong Zhang 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage); 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC); 2052401956bSBarry Smith 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal); 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal); 208014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal); 209421c37bdSBarry Smith 21019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType); 211014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool ); 212014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool ); 213014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC); 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC); 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool ); 216f5a88c2aSLois Curfman McInnes 2172591b318SToby Isaac PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC,PetscInt*); 218014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt); 219014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt); 220b35a507dSBarry Smith 221014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]); 222014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]); 223014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt); 224d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool); 225d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*); 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool ); 227f746d493SDmitry Karpeev 2283d957683SBarry Smith /*E 2293d957683SBarry Smith PCASMType - Type of additive Schwarz method to use 2303d957683SBarry Smith 231feb221f8SDmitry Karpeev $ PC_ASM_BASIC - Symmetric version where residuals from the ghost points are used 232feb221f8SDmitry Karpeev $ and computed values in ghost regions are added together. 233feb221f8SDmitry Karpeev $ Classical standard additive Schwarz. 234feb221f8SDmitry Karpeev $ PC_ASM_RESTRICT - Residuals from ghost points are used but computed values in ghost 235feb221f8SDmitry Karpeev $ region are discarded. 236feb221f8SDmitry Karpeev $ Default. 237feb221f8SDmitry Karpeev $ PC_ASM_INTERPOLATE - Residuals from ghost points are not used, computed values in ghost 238feb221f8SDmitry Karpeev $ region are added back in. 239feb221f8SDmitry Karpeev $ PC_ASM_NONE - Residuals from ghost points are not used, computed ghost values are 240feb221f8SDmitry Karpeev $ discarded. 241feb221f8SDmitry Karpeev $ Not very good. 2423d957683SBarry Smith 2433d957683SBarry Smith Level: beginner 2443d957683SBarry Smith 2453d957683SBarry Smith .seealso: PCASMSetType() 2463d957683SBarry Smith E*/ 247d252947aSBarry Smith typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 2486a6fc655SJed Brown PETSC_EXTERN const char *const PCASMTypes[]; 2499dcbbd2bSBarry Smith 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType); 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]); 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]); 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]); 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]); 256981c4779SBarry Smith 2573d957683SBarry Smith /*E 2586a4f0f73SDmitry Karpeev PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain). 259f746d493SDmitry Karpeev 2606a4f0f73SDmitry Karpeev Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 2616a4f0f73SDmitry Karpeev domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 2626a4f0f73SDmitry Karpeev a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 2636a4f0f73SDmitry Karpeev (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 2646a4f0f73SDmitry Karpeev 2656a4f0f73SDmitry Karpeev $ PC_GASM_BASIC - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 2666a4f0f73SDmitry Karpeev $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 2676a4f0f73SDmitry Karpeev $ from neighboring subdomains. 268feb221f8SDmitry Karpeev $ Classical standard additive Schwarz. 2696a4f0f73SDmitry Karpeev $ PC_GASM_RESTRICT - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 2706a4f0f73SDmitry Karpeev $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 2716a4f0f73SDmitry Karpeev $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 2726a4f0f73SDmitry Karpeev $ assumption). 273feb221f8SDmitry Karpeev $ Default. 2746a4f0f73SDmitry Karpeev $ PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 2756a4f0f73SDmitry Karpeev $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 2766a4f0f73SDmitry Karpeev $ from neighboring subdomains. 2776a4f0f73SDmitry Karpeev $ 2786a4f0f73SDmitry Karpeev $ PC_GASM_NONE - Residuals and corrections are zeroed out outside the local subdomains. 279feb221f8SDmitry Karpeev $ Not very good. 280f746d493SDmitry Karpeev 281f746d493SDmitry Karpeev Level: beginner 282f746d493SDmitry Karpeev 283f746d493SDmitry Karpeev .seealso: PCGASMSetType() 284f746d493SDmitry Karpeev E*/ 285f746d493SDmitry Karpeev typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType; 2866a6fc655SJed Brown PETSC_EXTERN const char *const PCGASMTypes[]; 287f746d493SDmitry Karpeev 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]); 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt); 291d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCGASMSetDMSubdomains(PC,PetscBool); 292d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCGASMGetDMSubdomains(PC,PetscBool*); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool ); 294f746d493SDmitry Karpeev 295014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType); 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]); 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]); 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 299014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]); 300014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]); 301f746d493SDmitry Karpeev 302f746d493SDmitry Karpeev /*E 3033d957683SBarry Smith PCCompositeType - Determines how two or more preconditioner are composed 3043d957683SBarry Smith 3053d957683SBarry Smith $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 3063d957683SBarry Smith $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 3073d957683SBarry Smith $ computed after the previous preconditioner application 308ccb205f8SBarry Smith $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 309ccb205f8SBarry Smith $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions) 3103d957683SBarry Smith $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 3113d957683SBarry Smith $ where first preconditioner is built from alpha I + S and second from 3123d957683SBarry Smith $ alpha I + R 3133d957683SBarry Smith 3143d957683SBarry Smith Level: beginner 3153d957683SBarry Smith 3163d957683SBarry Smith .seealso: PCCompositeSetType() 3173d957683SBarry Smith E*/ 3183b224e63SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType; 3196a6fc655SJed Brown PETSC_EXTERN const char *const PCCompositeTypes[]; 3209dcbbd2bSBarry Smith 321014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType); 32219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType); 323014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *); 324014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar); 325981c4779SBarry Smith 326014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt); 327014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter); 328014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*); 329da3a660dSBarry Smith 330014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double); 331014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt); 332014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt); 333014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt); 334014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt); 335014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt); 336014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt); 337014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt); 3383304466cSBarry Smith 339014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]); 340014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]); 341014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]); 342014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]); 3433304466cSBarry Smith 344014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*); 345014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*); 346014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType); 347014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt); 348014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS); 349014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*); 3504ab8060aSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool); 3514ab8060aSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*); 352084e4875SJed Brown 353084e4875SJed Brown /*E 354084e4875SJed Brown PCFieldSplitSchurPreType - Determines how to precondition Schur complement 355084e4875SJed Brown 356084e4875SJed Brown Level: intermediate 357084e4875SJed Brown 358084e4875SJed Brown .seealso: PCFieldSplitSchurPrecondition() 359084e4875SJed Brown E*/ 360*2e71c61dSMatthew G. Knepley typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_SELFP,PC_FIELDSPLIT_SCHUR_PRE_A11,PC_FIELDSPLIT_SCHUR_PRE_USER,PC_FIELDSPLIT_SCHUR_PRE_FULL} PCFieldSplitSchurPreType; 361014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[]; 362084e4875SJed Brown 363ab1df9f5SJed Brown /*E 364c9c6ffaaSJed Brown PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 365ab1df9f5SJed Brown 366ab1df9f5SJed Brown Level: intermediate 367ab1df9f5SJed Brown 368c9c6ffaaSJed Brown .seealso: PCFieldSplitSetSchurFactType() 369ab1df9f5SJed Brown E*/ 370ab1df9f5SJed Brown typedef enum { 371c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_DIAG, 372c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_LOWER, 373c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_UPPER, 374c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_FULL 375c9c6ffaaSJed Brown } PCFieldSplitSchurFactType; 376014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[]; 377ab1df9f5SJed Brown 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat); 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*); 381470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC,Mat *S); 382470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC,Mat *S); 3833d30b1ffSBarry Smith 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat); 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat); 3862a6744ebSBarry Smith 387014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*); 388014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *); 3892102ba4dSSatish Balay 390014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]); 39167fac13cSBarry Smith 392014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM); 393014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*); 3946c699258SBarry Smith 395014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*); 396014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*); 3971b2093e4SBarry Smith 398014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal); 399014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt); 400014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool); 4016c699258SBarry Smith 402014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal); 403014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool); 404014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt); 405014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt); 40637f80224SJose E. Roman /*E 40737f80224SJose E. Roman PCPARMSGlobalType - Determines the global preconditioner method in PARMS 40837f80224SJose E. Roman 40937f80224SJose E. Roman Level: intermediate 41037f80224SJose E. Roman 41137f80224SJose E. Roman .seealso: PCPARMSSetGlobal() 41237f80224SJose E. Roman E*/ 41337f80224SJose E. Roman typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType; 4146a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSGlobalTypes[]; 41537f80224SJose E. Roman /*E 41637f80224SJose E. Roman PCPARMSLocalType - Determines the local preconditioner method in PARMS 41737f80224SJose E. Roman 41837f80224SJose E. Roman Level: intermediate 41937f80224SJose E. Roman 42037f80224SJose E. Roman .seealso: PCPARMSSetLocal() 42137f80224SJose E. Roman E*/ 42237f80224SJose E. Roman typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType; 4236a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSLocalTypes[]; 42437f80224SJose E. Roman 425014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type); 426014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type); 427014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits); 428014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart); 429014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym); 430014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2); 43137f80224SJose E. Roman 4325154939eSJed Brown /*E 4335154939eSJed Brown PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method 4345154939eSJed Brown 4355154939eSJed Brown Level: intermediate 4365154939eSJed Brown 4375154939eSJed Brown .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseProl() 4385154939eSJed Brown E*/ 439a782d3a5SMark Adams typedef const char *PCGAMGType; 440a782d3a5SMark Adams #define PCGAMGAGG "agg" 441a782d3a5SMark Adams #define PCGAMGGEO "geo" 4428e6d0c30SPeter Brune #define PCGAMGCLASSICAL "classical" 443014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt); 444014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool); 445014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool); 446014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt); 447014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal); 448014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt); 449014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt); 45019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType ); 451014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n); 452014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n); 453014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool); 454fc897844SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool); 4553e3471ccSMark Adams PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void); 4563e3471ccSMark Adams PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void); 45756de70efSSatish Balay 4588eab0cc1SPeter Brune typedef const char *PCGAMGClassicalType; 4598eab0cc1SPeter Brune #define PCGAMGCLASSICALDIRECT "direct" 4608eab0cc1SPeter Brune #define PCGAMGCLASSICALSTANDARD "standard" 4618eab0cc1SPeter Brune PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC,PCGAMGClassicalType); 4628eab0cc1SPeter Brune 4630c7d97c5SJed Brown #if defined(PETSC_HAVE_PCBDDC) 464674ae819SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC,IS); 4654fad6a16SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt); 4662b510759SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC,PetscInt); 4670bdf917eSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace); 468014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS); 469014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*); 470014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS); 471014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*); 472014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]); 4731a83f524SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode); 4743425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*); 4753425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec); 4763425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec); 4770c7d97c5SJed Brown #endif 4780c7d97c5SJed Brown 479b965d45eSStefano Zampini PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool); 480014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar); 48146caae47SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec); 482d3b1e0e7SMatthew Knepley 483b4876bcbSBarry Smith /*E 484b4876bcbSBarry Smith PCMGType - Determines the type of multigrid method that is run. 485b4876bcbSBarry Smith 486b4876bcbSBarry Smith Level: beginner 487b4876bcbSBarry Smith 488b4876bcbSBarry Smith Values: 489b4876bcbSBarry Smith + PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles() 490b4876bcbSBarry Smith . PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are 491b4876bcbSBarry Smith smoothed before updating the residual. This only uses the 492b4876bcbSBarry Smith down smoother, in the preconditioner the upper smoother is ignored 493b4876bcbSBarry Smith . PC_MG_FULL - same as multiplicative except one also performs grid sequencing, 494b4876bcbSBarry Smith that is starts on the coarsest grid, performs a cycle, interpolates 495b4876bcbSBarry 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 496b4876bcbSBarry Smith algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 497b4876bcbSBarry Smith calls the V-cycle only on the coarser level and has a post-smoother instead. 498b4876bcbSBarry Smith - PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level 499b4876bcbSBarry Smith from a finer 500b4876bcbSBarry Smith 501b4876bcbSBarry Smith .seealso: PCMGSetType() 502b4876bcbSBarry Smith 503b4876bcbSBarry Smith E*/ 504b4876bcbSBarry Smith typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType; 505b4876bcbSBarry Smith PETSC_EXTERN const char *const PCMGTypes[]; 506b4876bcbSBarry Smith #define PC_MG_CASCADE PC_MG_KASKADE; 507b4876bcbSBarry Smith 508b4876bcbSBarry Smith /*E 509b4876bcbSBarry Smith PCMGCycleType - Use V-cycle or W-cycle 510b4876bcbSBarry Smith 511b4876bcbSBarry Smith Level: beginner 512b4876bcbSBarry Smith 513b4876bcbSBarry Smith Values: 514b4876bcbSBarry Smith + PC_MG_V_CYCLE 515b4876bcbSBarry Smith - PC_MG_W_CYCLE 516b4876bcbSBarry Smith 517b4876bcbSBarry Smith .seealso: PCMGSetCycleType() 518b4876bcbSBarry Smith 519b4876bcbSBarry Smith E*/ 520b4876bcbSBarry Smith typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType; 521b4876bcbSBarry Smith PETSC_EXTERN const char *const PCMGCycleTypes[]; 522b4876bcbSBarry Smith 523b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType); 524b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*); 525b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*); 526b4876bcbSBarry Smith 527b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt); 528b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt); 529b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType); 530b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType); 531b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt); 532b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt); 533b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool); 534b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*); 535b4876bcbSBarry Smith 536b4876bcbSBarry Smith 537b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec); 538b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec); 539b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec); 540b4876bcbSBarry Smith 541b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat); 542b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*); 543b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat); 544b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*); 545b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec); 546b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*); 547b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat); 54854b2cd4bSJed Brown PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat,Vec,Vec,Vec); 549b4876bcbSBarry Smith 550b4876bcbSBarry Smith /*E 551b4876bcbSBarry Smith PCExoticType - Face based or wirebasket based coarse grid space 552b4876bcbSBarry Smith 553b4876bcbSBarry Smith Level: beginner 554b4876bcbSBarry Smith 555b4876bcbSBarry Smith .seealso: PCExoticSetType(), PCEXOTIC 556b4876bcbSBarry Smith E*/ 557b4876bcbSBarry Smith typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType; 558b4876bcbSBarry Smith PETSC_EXTERN const char *const PCExoticTypes[]; 559b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType); 560b4876bcbSBarry Smith 561123ea438SMatthew Knepley #endif /* __PETSCPC_H */ 562