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" 6724c02a0fSBarry Smith #define PCCP "cp" 68628b8437SMatthew Knepley #define PCBFBT "bfbt" 69519d70e2SJed Brown #define PCLSC "lsc" 701d6018f0SLisandro Dalcin #define PCPYTHON "python" 71f91d8e95SBarry Smith #define PCPFMG "pfmg" 72d851a50bSGlenn Hammond #define PCSYSPFMG "syspfmg" 73df826632SBarry Smith #define PCREDISTRIBUTE "redistribute" 74ac793be5SBarry Smith #define PCSVD "svd" 75ac793be5SBarry Smith #define PCGAMG "gamg" 76ac793be5SBarry Smith #define PCSACUSP "sacusp" /* these four run on NVIDIA GPUs using CUSP */ 77bfc29b71SVictor Minden #define PCSACUSPPOLY "sacusppoly" 788154be41SBarry Smith #define PCBICGSTABCUSP "bicgstabcusp" 7904b59e88SVictor Minden #define PCAINVCUSP "ainvcusp" 800c7d97c5SJed Brown #define PCBDDC "bddc" 813f93e5bdSPeter Brune #define PCKACZMARZ "kaczmarz" 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); 10923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode PCSetReusePreconditioner(PC,PetscBool); 1104d0a8057SBarry Smith 11155849f57SBarry Smith #define PC_FILE_CLASSID 1211222 11255849f57SBarry Smith 1134d0a8057SBarry Smith /*E 1144d0a8057SBarry Smith PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates 1154d0a8057SBarry Smith 1164d0a8057SBarry Smith Level: advanced 1174d0a8057SBarry Smith 1184d0a8057SBarry Smith Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h 1194d0a8057SBarry Smith 1204d0a8057SBarry Smith .seealso: PCApplyRichardson() 1214d0a8057SBarry Smith E*/ 1224d0a8057SBarry Smith typedef enum { 1234d0a8057SBarry Smith PCRICHARDSON_CONVERGED_RTOL = 2, 1244d0a8057SBarry Smith PCRICHARDSON_CONVERGED_ATOL = 3, 1254d0a8057SBarry Smith PCRICHARDSON_CONVERGED_ITS = 4, 1264d0a8057SBarry Smith PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason; 1274d0a8057SBarry Smith 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*); 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *); 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool ); 13149517cdeSBarry Smith PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool); 13249517cdeSBarry Smith PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*); 13384cb2905SBarry Smith 134607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PCRegisterAll(void); 135014dd563SJed Brown PETSC_EXTERN PetscBool PCRegisterAllCalled; 13684cb2905SBarry Smith 137bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode PCRegister(const char[],PetscErrorCode(*)(PC)); 13830de9b25SBarry Smith 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCReset(PC); 140014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDestroy(PC*); 141014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC); 14219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*); 14314c91fddSBarry Smith 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*); 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*); 146014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*); 1475b116368SBarry Smith 14823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat); 14923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*); 150014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *); 1514b0e389bSBarry Smith 152014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer); 15355849f57SBarry Smith PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer); 154ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PCViewFromOptions(PC A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 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 173*baa89ecbSBarry Smith /*E 174*baa89ecbSBarry Smith PCJacobiType - What elements are used to form the Jacobi preconditioner 175*baa89ecbSBarry Smith 176*baa89ecbSBarry Smith Level: intermediate 177*baa89ecbSBarry Smith 178*baa89ecbSBarry Smith .seealso: 179*baa89ecbSBarry Smith E*/ 180*baa89ecbSBarry Smith typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType; 181*baa89ecbSBarry Smith PETSC_EXTERN const char *const PCJacobiTypes[]; 182*baa89ecbSBarry Smith 183*baa89ecbSBarry Smith PETSC_EXTERN PetscErrorCode PCJacobiSetType(PC,PCJacobiType); 184*baa89ecbSBarry Smith PETSC_EXTERN PetscErrorCode PCJacobiGetType(PC,PCJacobiType*); 185*baa89ecbSBarry Smith PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC,PetscBool); 186*baa89ecbSBarry Smith PETSC_EXTERN PetscErrorCode PCJacobiGetUseAbs(PC,PetscBool*); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType); 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt); 190d03aef70SBarry Smith 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal); 192014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC); 193421c37bdSBarry Smith 194014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]); 195014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]); 1961eb62cbbSBarry Smith 197014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec)); 198014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)); 199014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec)); 200014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC)); 201014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*)); 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer)); 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC)); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*); 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]); 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]); 208aabeff55SBarry Smith 209014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal); 210d90ac83dSHong Zhang 211014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType); 212014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal); 213d90ac83dSHong Zhang 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage); 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*); 216014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC); 2172401956bSBarry Smith 218014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal); 219014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal); 220014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal); 221421c37bdSBarry Smith 22219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType); 223014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool ); 224014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool ); 225014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC); 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC); 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool ); 228f5a88c2aSLois Curfman McInnes 2292591b318SToby Isaac PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC,PetscInt*); 230014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt); 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt); 232b35a507dSBarry Smith 233014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]); 234014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]); 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt); 236d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool); 237d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*); 238014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool ); 239f746d493SDmitry Karpeev 2403d957683SBarry Smith /*E 2413d957683SBarry Smith PCASMType - Type of additive Schwarz method to use 2423d957683SBarry Smith 243feb221f8SDmitry Karpeev $ PC_ASM_BASIC - Symmetric version where residuals from the ghost points are used 244feb221f8SDmitry Karpeev $ and computed values in ghost regions are added together. 245feb221f8SDmitry Karpeev $ Classical standard additive Schwarz. 246feb221f8SDmitry Karpeev $ PC_ASM_RESTRICT - Residuals from ghost points are used but computed values in ghost 247feb221f8SDmitry Karpeev $ region are discarded. 248feb221f8SDmitry Karpeev $ Default. 249feb221f8SDmitry Karpeev $ PC_ASM_INTERPOLATE - Residuals from ghost points are not used, computed values in ghost 250feb221f8SDmitry Karpeev $ region are added back in. 251feb221f8SDmitry Karpeev $ PC_ASM_NONE - Residuals from ghost points are not used, computed ghost values are 252feb221f8SDmitry Karpeev $ discarded. 253feb221f8SDmitry Karpeev $ Not very good. 2543d957683SBarry Smith 2553d957683SBarry Smith Level: beginner 2563d957683SBarry Smith 2573d957683SBarry Smith .seealso: PCASMSetType() 2583d957683SBarry Smith E*/ 259d252947aSBarry Smith typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 2606a6fc655SJed Brown PETSC_EXTERN const char *const PCASMTypes[]; 2619dcbbd2bSBarry Smith 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 266014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]); 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]); 268981c4779SBarry Smith 2693d957683SBarry Smith /*E 2706a4f0f73SDmitry Karpeev PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain). 271f746d493SDmitry Karpeev 2726a4f0f73SDmitry Karpeev Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 2736a4f0f73SDmitry Karpeev domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 2746a4f0f73SDmitry Karpeev a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 2756a4f0f73SDmitry Karpeev (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 2766a4f0f73SDmitry Karpeev 2776a4f0f73SDmitry Karpeev $ PC_GASM_BASIC - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 2786a4f0f73SDmitry Karpeev $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 2796a4f0f73SDmitry Karpeev $ from neighboring subdomains. 280feb221f8SDmitry Karpeev $ Classical standard additive Schwarz. 2816a4f0f73SDmitry Karpeev $ PC_GASM_RESTRICT - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 2826a4f0f73SDmitry Karpeev $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 2836a4f0f73SDmitry Karpeev $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 2846a4f0f73SDmitry Karpeev $ assumption). 285feb221f8SDmitry Karpeev $ Default. 2866a4f0f73SDmitry Karpeev $ PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 2876a4f0f73SDmitry Karpeev $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 2886a4f0f73SDmitry Karpeev $ from neighboring subdomains. 2896a4f0f73SDmitry Karpeev $ 2906a4f0f73SDmitry Karpeev $ PC_GASM_NONE - Residuals and corrections are zeroed out outside the local subdomains. 291feb221f8SDmitry Karpeev $ Not very good. 292f746d493SDmitry Karpeev 293f746d493SDmitry Karpeev Level: beginner 294f746d493SDmitry Karpeev 295f746d493SDmitry Karpeev .seealso: PCGASMSetType() 296f746d493SDmitry Karpeev E*/ 297f746d493SDmitry Karpeev typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType; 2986a6fc655SJed Brown PETSC_EXTERN const char *const PCGASMTypes[]; 299f746d493SDmitry Karpeev 300014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]); 301014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool); 302014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt); 303d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCGASMSetDMSubdomains(PC,PetscBool); 304d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCGASMGetDMSubdomains(PC,PetscBool*); 305014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool ); 306f746d493SDmitry Karpeev 307014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType); 308014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]); 309014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]); 310014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 311014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]); 312014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]); 313f746d493SDmitry Karpeev 314f746d493SDmitry Karpeev /*E 3153d957683SBarry Smith PCCompositeType - Determines how two or more preconditioner are composed 3163d957683SBarry Smith 3173d957683SBarry Smith $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 3183d957683SBarry Smith $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 3193d957683SBarry Smith $ computed after the previous preconditioner application 320ccb205f8SBarry Smith $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 321ccb205f8SBarry Smith $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions) 3223d957683SBarry Smith $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 3233d957683SBarry Smith $ where first preconditioner is built from alpha I + S and second from 3243d957683SBarry Smith $ alpha I + R 3253d957683SBarry Smith 3263d957683SBarry Smith Level: beginner 3273d957683SBarry Smith 3283d957683SBarry Smith .seealso: PCCompositeSetType() 3293d957683SBarry Smith E*/ 3303b224e63SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType; 3316a6fc655SJed Brown PETSC_EXTERN const char *const PCCompositeTypes[]; 3329dcbbd2bSBarry Smith 333014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType); 33419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType); 335014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *); 336014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar); 337981c4779SBarry Smith 338014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt); 339014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter); 340014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*); 341da3a660dSBarry Smith 342014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double); 343014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt); 344014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt); 345014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt); 346014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt); 347014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt); 348014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt); 349014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt); 3503304466cSBarry Smith 351014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]); 352014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]); 353014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]); 354014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]); 3553304466cSBarry Smith 356014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*); 357014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*); 358014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType); 359014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt); 360014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS); 361014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*); 3624ab8060aSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool); 3634ab8060aSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*); 364c84da90fSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC,PetscBool); 365c84da90fSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC,PetscBool*); 366c84da90fSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC,PetscBool); 367c84da90fSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC,PetscBool*); 368c84da90fSDmitry Karpeev 369084e4875SJed Brown 370084e4875SJed Brown /*E 371084e4875SJed Brown PCFieldSplitSchurPreType - Determines how to precondition Schur complement 372084e4875SJed Brown 373084e4875SJed Brown Level: intermediate 374084e4875SJed Brown 37529f8a81cSJed Brown .seealso: PCFieldSplitSetSchurPre() 376084e4875SJed Brown E*/ 3772e71c61dSMatthew 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; 378014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[]; 379084e4875SJed Brown 380ab1df9f5SJed Brown /*E 381c9c6ffaaSJed Brown PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 382ab1df9f5SJed Brown 383ab1df9f5SJed Brown Level: intermediate 384ab1df9f5SJed Brown 385c9c6ffaaSJed Brown .seealso: PCFieldSplitSetSchurFactType() 386ab1df9f5SJed Brown E*/ 387ab1df9f5SJed Brown typedef enum { 388c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_DIAG, 389c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_LOWER, 390c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_UPPER, 391c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_FULL 392c9c6ffaaSJed Brown } PCFieldSplitSchurFactType; 393014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[]; 394ab1df9f5SJed Brown 39529f8a81cSJed Brown PETSC_EXTERN PETSC_DEPRECATED("Use PCFieldSplitSetSchurPre") PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat); 39629f8a81cSJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC,PCFieldSplitSchurPreType,Mat); 39737a82bf0SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC,PCFieldSplitSchurPreType*,Mat*); 398014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType); 399014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*); 400470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC,Mat *S); 401470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC,Mat *S); 4023d30b1ffSBarry Smith 403014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat); 404014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat); 4052a6744ebSBarry Smith 406014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*); 4072102ba4dSSatish Balay 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]); 40967fac13cSBarry Smith 410014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM); 411014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*); 4126c699258SBarry Smith 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*); 414014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*); 4151b2093e4SBarry Smith 416014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal); 417014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt); 418014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool); 4196c699258SBarry Smith 420014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal); 421014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool); 422014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt); 423014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt); 42437f80224SJose E. Roman /*E 42537f80224SJose E. Roman PCPARMSGlobalType - Determines the global preconditioner method in PARMS 42637f80224SJose E. Roman 42737f80224SJose E. Roman Level: intermediate 42837f80224SJose E. Roman 42937f80224SJose E. Roman .seealso: PCPARMSSetGlobal() 43037f80224SJose E. Roman E*/ 43137f80224SJose E. Roman typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType; 4326a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSGlobalTypes[]; 43337f80224SJose E. Roman /*E 43437f80224SJose E. Roman PCPARMSLocalType - Determines the local preconditioner method in PARMS 43537f80224SJose E. Roman 43637f80224SJose E. Roman Level: intermediate 43737f80224SJose E. Roman 43837f80224SJose E. Roman .seealso: PCPARMSSetLocal() 43937f80224SJose E. Roman E*/ 44037f80224SJose E. Roman typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType; 4416a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSLocalTypes[]; 44237f80224SJose E. Roman 443*baa89ecbSBarry Smith PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC,PCPARMSGlobalType); 444*baa89ecbSBarry Smith PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC,PCPARMSLocalType); 445*baa89ecbSBarry Smith PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC,PetscReal,PetscInt); 446*baa89ecbSBarry Smith PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC,PetscInt); 447*baa89ecbSBarry Smith PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC,PetscBool); 448*baa89ecbSBarry Smith PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC,PetscInt,PetscInt,PetscInt); 44937f80224SJose E. Roman 4505154939eSJed Brown /*E 4515154939eSJed Brown PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method 4525154939eSJed Brown 4535154939eSJed Brown Level: intermediate 4545154939eSJed Brown 4555154939eSJed Brown .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseProl() 4565154939eSJed Brown E*/ 457a782d3a5SMark Adams typedef const char *PCGAMGType; 458a782d3a5SMark Adams #define PCGAMGAGG "agg" 459a782d3a5SMark Adams #define PCGAMGGEO "geo" 4608e6d0c30SPeter Brune #define PCGAMGCLASSICAL "classical" 461014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt); 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool); 463014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool); 464014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt); 465014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal); 466014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt); 467014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt); 46819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType); 469*baa89ecbSBarry Smith PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC,PetscInt); 470*baa89ecbSBarry Smith PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC,PetscBool); 471014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool); 472fc897844SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool); 4733e3471ccSMark Adams PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void); 4743e3471ccSMark Adams PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void); 47556de70efSSatish Balay 4768eab0cc1SPeter Brune typedef const char *PCGAMGClassicalType; 4778eab0cc1SPeter Brune #define PCGAMGCLASSICALDIRECT "direct" 4788eab0cc1SPeter Brune #define PCGAMGCLASSICALSTANDARD "standard" 4798eab0cc1SPeter Brune PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC,PCGAMGClassicalType); 4808eab0cc1SPeter Brune 4810c7d97c5SJed Brown #if defined(PETSC_HAVE_PCBDDC) 482906d46d4SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisMat(PC,Mat); 483674ae819SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC,IS); 4844fad6a16SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt); 4852b510759SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC,PetscInt); 4860bdf917eSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace); 487014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS); 48882ba6b80SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC,IS); 489014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*); 49082ba6b80SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC,IS*); 491014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS); 49282ba6b80SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC,IS); 493014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*); 49482ba6b80SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC,IS*); 495014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]); 49663602bcaSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC,PetscInt,IS[]); 4971a83f524SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode); 4983425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*); 4993425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec); 5003425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec); 5010c7d97c5SJed Brown #endif 5020c7d97c5SJed Brown 503b965d45eSStefano Zampini PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool); 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar); 50546caae47SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec); 506d3b1e0e7SMatthew Knepley 507b4876bcbSBarry Smith /*E 508b4876bcbSBarry Smith PCMGType - Determines the type of multigrid method that is run. 509b4876bcbSBarry Smith 510b4876bcbSBarry Smith Level: beginner 511b4876bcbSBarry Smith 512b4876bcbSBarry Smith Values: 513b4876bcbSBarry Smith + PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles() 514b4876bcbSBarry Smith . PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are 515b4876bcbSBarry Smith smoothed before updating the residual. This only uses the 516b4876bcbSBarry Smith down smoother, in the preconditioner the upper smoother is ignored 517b4876bcbSBarry Smith . PC_MG_FULL - same as multiplicative except one also performs grid sequencing, 518b4876bcbSBarry Smith that is starts on the coarsest grid, performs a cycle, interpolates 519b4876bcbSBarry 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 520b4876bcbSBarry Smith algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 521b4876bcbSBarry Smith calls the V-cycle only on the coarser level and has a post-smoother instead. 522b4876bcbSBarry Smith - PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level 523b4876bcbSBarry Smith from a finer 524b4876bcbSBarry Smith 525b4876bcbSBarry Smith .seealso: PCMGSetType() 526b4876bcbSBarry Smith 527b4876bcbSBarry Smith E*/ 528b4876bcbSBarry Smith typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType; 529b4876bcbSBarry Smith PETSC_EXTERN const char *const PCMGTypes[]; 530b4876bcbSBarry Smith #define PC_MG_CASCADE PC_MG_KASKADE; 531b4876bcbSBarry Smith 532b4876bcbSBarry Smith /*E 533b4876bcbSBarry Smith PCMGCycleType - Use V-cycle or W-cycle 534b4876bcbSBarry Smith 535b4876bcbSBarry Smith Level: beginner 536b4876bcbSBarry Smith 537b4876bcbSBarry Smith Values: 538b4876bcbSBarry Smith + PC_MG_V_CYCLE 539b4876bcbSBarry Smith - PC_MG_W_CYCLE 540b4876bcbSBarry Smith 541b4876bcbSBarry Smith .seealso: PCMGSetCycleType() 542b4876bcbSBarry Smith 543b4876bcbSBarry Smith E*/ 544b4876bcbSBarry Smith typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType; 545b4876bcbSBarry Smith PETSC_EXTERN const char *const PCMGCycleTypes[]; 546b4876bcbSBarry Smith 547b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType); 548b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*); 549b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*); 550b4876bcbSBarry Smith 551b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt); 552b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt); 553b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType); 554b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType); 555b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt); 556b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt); 557b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool); 558b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*); 559b4876bcbSBarry Smith 560b4876bcbSBarry Smith 561b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec); 562b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec); 563b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec); 564b4876bcbSBarry Smith 565b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat); 566b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*); 567b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat); 568b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*); 569b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec); 570b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*); 571b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat); 57254b2cd4bSJed Brown PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat,Vec,Vec,Vec); 573b4876bcbSBarry Smith 574b4876bcbSBarry Smith /*E 575b4876bcbSBarry Smith PCExoticType - Face based or wirebasket based coarse grid space 576b4876bcbSBarry Smith 577b4876bcbSBarry Smith Level: beginner 578b4876bcbSBarry Smith 579b4876bcbSBarry Smith .seealso: PCExoticSetType(), PCEXOTIC 580b4876bcbSBarry Smith E*/ 581b4876bcbSBarry Smith typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType; 582b4876bcbSBarry Smith PETSC_EXTERN const char *const PCExoticTypes[]; 583b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType); 584b4876bcbSBarry Smith 585123ea438SMatthew Knepley #endif /* __PETSCPC_H */ 586