xref: /petsc/include/petscpc.h (revision 906d46d482912ff65fe95caf6add84f15a3f7e20)
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 
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*);
353c84da90fSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC,PetscBool);
354c84da90fSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC,PetscBool*);
355c84da90fSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC,PetscBool);
356c84da90fSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC,PetscBool*);
357c84da90fSDmitry Karpeev 
358084e4875SJed Brown 
359084e4875SJed Brown /*E
360084e4875SJed Brown     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
361084e4875SJed Brown 
362084e4875SJed Brown     Level: intermediate
363084e4875SJed Brown 
36429f8a81cSJed Brown .seealso: PCFieldSplitSetSchurPre()
365084e4875SJed Brown E*/
3662e71c61dSMatthew 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;
367014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
368084e4875SJed Brown 
369ab1df9f5SJed Brown /*E
370c9c6ffaaSJed Brown     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
371ab1df9f5SJed Brown 
372ab1df9f5SJed Brown     Level: intermediate
373ab1df9f5SJed Brown 
374c9c6ffaaSJed Brown .seealso: PCFieldSplitSetSchurFactType()
375ab1df9f5SJed Brown E*/
376ab1df9f5SJed Brown typedef enum {
377c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
378c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
379c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
380c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_FULL
381c9c6ffaaSJed Brown } PCFieldSplitSchurFactType;
382014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
383ab1df9f5SJed Brown 
38429f8a81cSJed Brown PETSC_EXTERN PETSC_DEPRECATED("Use PCFieldSplitSetSchurPre") PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
38529f8a81cSJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC,PCFieldSplitSchurPreType,Mat);
38637a82bf0SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC,PCFieldSplitSchurPreType*,Mat*);
387014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
388014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
389470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC,Mat *S);
390470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC,Mat *S);
3913d30b1ffSBarry Smith 
392014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
393014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
3942a6744ebSBarry Smith 
395014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
396014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);
3972102ba4dSSatish Balay 
398014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);
39967fac13cSBarry Smith 
400014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
401014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);
4026c699258SBarry Smith 
403014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);
4051b2093e4SBarry Smith 
406014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
407014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
408014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
4096c699258SBarry Smith 
410014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
411014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
412014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
413014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
41437f80224SJose E. Roman /*E
41537f80224SJose E. Roman     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
41637f80224SJose E. Roman 
41737f80224SJose E. Roman     Level: intermediate
41837f80224SJose E. Roman 
41937f80224SJose E. Roman .seealso: PCPARMSSetGlobal()
42037f80224SJose E. Roman E*/
42137f80224SJose E. Roman typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
4226a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
42337f80224SJose E. Roman /*E
42437f80224SJose E. Roman     PCPARMSLocalType - Determines the local preconditioner method in PARMS
42537f80224SJose E. Roman 
42637f80224SJose E. Roman     Level: intermediate
42737f80224SJose E. Roman 
42837f80224SJose E. Roman .seealso: PCPARMSSetLocal()
42937f80224SJose E. Roman E*/
43037f80224SJose E. Roman typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
4316a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSLocalTypes[];
43237f80224SJose E. Roman 
433014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type);
434014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type);
435014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits);
436014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart);
437014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym);
438014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2);
43937f80224SJose E. Roman 
4405154939eSJed Brown /*E
4415154939eSJed Brown     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method
4425154939eSJed Brown 
4435154939eSJed Brown     Level: intermediate
4445154939eSJed Brown 
4455154939eSJed Brown .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseProl()
4465154939eSJed Brown E*/
447a782d3a5SMark Adams typedef const char *PCGAMGType;
448a782d3a5SMark Adams #define PCGAMGAGG         "agg"
449a782d3a5SMark Adams #define PCGAMGGEO         "geo"
4508e6d0c30SPeter Brune #define PCGAMGCLASSICAL   "classical"
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
453014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
454014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
455014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
456014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
457014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
45819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType );
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n);
460014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n);
461014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool);
462fc897844SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool);
4633e3471ccSMark Adams PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void);
4643e3471ccSMark Adams PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void);
46556de70efSSatish Balay 
4668eab0cc1SPeter Brune typedef const char *PCGAMGClassicalType;
4678eab0cc1SPeter Brune #define PCGAMGCLASSICALDIRECT   "direct"
4688eab0cc1SPeter Brune #define PCGAMGCLASSICALSTANDARD "standard"
4698eab0cc1SPeter Brune PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC,PCGAMGClassicalType);
4708eab0cc1SPeter Brune 
4710c7d97c5SJed Brown #if defined(PETSC_HAVE_PCBDDC)
472*906d46d4SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisMat(PC,Mat);
473674ae819SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC,IS);
4744fad6a16SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
4752b510759SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC,PetscInt);
4760bdf917eSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace);
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
47882ba6b80SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC,IS);
479014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
48082ba6b80SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC,IS*);
481014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
48282ba6b80SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC,IS);
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
48482ba6b80SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC,IS*);
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
48663602bcaSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC,PetscInt,IS[]);
4871a83f524SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode);
4883425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*);
4893425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec);
4903425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec);
4910c7d97c5SJed Brown #endif
4920c7d97c5SJed Brown 
493b965d45eSStefano Zampini PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool);
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
49546caae47SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec);
496d3b1e0e7SMatthew Knepley 
497b4876bcbSBarry Smith /*E
498b4876bcbSBarry Smith     PCMGType - Determines the type of multigrid method that is run.
499b4876bcbSBarry Smith 
500b4876bcbSBarry Smith    Level: beginner
501b4876bcbSBarry Smith 
502b4876bcbSBarry Smith    Values:
503b4876bcbSBarry Smith +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles()
504b4876bcbSBarry Smith .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
505b4876bcbSBarry Smith                 smoothed before updating the residual. This only uses the
506b4876bcbSBarry Smith                 down smoother, in the preconditioner the upper smoother is ignored
507b4876bcbSBarry Smith .  PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
508b4876bcbSBarry Smith             that is starts on the coarsest grid, performs a cycle, interpolates
509b4876bcbSBarry 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
510b4876bcbSBarry Smith             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
511b4876bcbSBarry Smith             calls the V-cycle only on the coarser level and has a post-smoother instead.
512b4876bcbSBarry Smith -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
513b4876bcbSBarry Smith                from a finer
514b4876bcbSBarry Smith 
515b4876bcbSBarry Smith .seealso: PCMGSetType()
516b4876bcbSBarry Smith 
517b4876bcbSBarry Smith E*/
518b4876bcbSBarry Smith typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
519b4876bcbSBarry Smith PETSC_EXTERN const char *const PCMGTypes[];
520b4876bcbSBarry Smith #define PC_MG_CASCADE PC_MG_KASKADE;
521b4876bcbSBarry Smith 
522b4876bcbSBarry Smith /*E
523b4876bcbSBarry Smith     PCMGCycleType - Use V-cycle or W-cycle
524b4876bcbSBarry Smith 
525b4876bcbSBarry Smith    Level: beginner
526b4876bcbSBarry Smith 
527b4876bcbSBarry Smith    Values:
528b4876bcbSBarry Smith +  PC_MG_V_CYCLE
529b4876bcbSBarry Smith -  PC_MG_W_CYCLE
530b4876bcbSBarry Smith 
531b4876bcbSBarry Smith .seealso: PCMGSetCycleType()
532b4876bcbSBarry Smith 
533b4876bcbSBarry Smith E*/
534b4876bcbSBarry Smith typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
535b4876bcbSBarry Smith PETSC_EXTERN const char *const PCMGCycleTypes[];
536b4876bcbSBarry Smith 
537b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType);
538b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*);
539b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*);
540b4876bcbSBarry Smith 
541b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt);
542b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt);
543b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType);
544b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType);
545b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt);
546b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt);
547b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool);
548b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*);
549b4876bcbSBarry Smith 
550b4876bcbSBarry Smith 
551b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec);
552b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec);
553b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec);
554b4876bcbSBarry Smith 
555b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat);
556b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*);
557b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat);
558b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*);
559b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec);
560b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*);
561b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
56254b2cd4bSJed Brown PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat,Vec,Vec,Vec);
563b4876bcbSBarry Smith 
564b4876bcbSBarry Smith /*E
565b4876bcbSBarry Smith     PCExoticType - Face based or wirebasket based coarse grid space
566b4876bcbSBarry Smith 
567b4876bcbSBarry Smith    Level: beginner
568b4876bcbSBarry Smith 
569b4876bcbSBarry Smith .seealso: PCExoticSetType(), PCEXOTIC
570b4876bcbSBarry Smith E*/
571b4876bcbSBarry Smith typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
572b4876bcbSBarry Smith PETSC_EXTERN const char *const PCExoticTypes[];
573b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);
574b4876bcbSBarry Smith 
575123ea438SMatthew Knepley #endif /* __PETSCPC_H */
576