xref: /petsc/include/petscpc.h (revision 014dd563d73e9fc78d056590fa6cf997782bf92d)
1d03aef70SBarry Smith /*
237f753daSBarry Smith       Preconditioner module.
3d03aef70SBarry Smith */
40a835dfdSSatish Balay #if !defined(__PETSCPC_H)
50a835dfdSSatish Balay #define __PETSCPC_H
6e1589f56SBarry Smith #include "petscdm.h"
7d03aef70SBarry Smith 
8*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCInitializePackage(const char[]);
91dbb0a54SBarry Smith 
10eec0b4cfSBarry Smith /*
11eec0b4cfSBarry Smith     PCList contains the list of preconditioners currently registered
12f1af5d2fSBarry Smith    These are added with the PCRegisterDynamic() macro
13eec0b4cfSBarry Smith */
14*014dd563SJed Brown PETSC_EXTERN PetscFList PCList;
1582bf6240SBarry Smith 
163d957683SBarry Smith /*S
173d957683SBarry Smith      PC - Abstract PETSc object that manages all preconditioners
183d957683SBarry Smith 
193d957683SBarry Smith    Level: beginner
203d957683SBarry Smith 
213d957683SBarry Smith   Concepts: preconditioners
223d957683SBarry Smith 
231a480d89SAdministrator .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
243d957683SBarry Smith S*/
253d957683SBarry Smith typedef struct _p_PC* PC;
263d957683SBarry Smith 
2776bdecfbSBarry Smith /*J
283d957683SBarry Smith     PCType - String with the name of a PETSc preconditioner method or the creation function
293d957683SBarry Smith        with an optional dynamic library name, for example
303d957683SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mypccreate()
313d957683SBarry Smith 
323d957683SBarry Smith    Level: beginner
333d957683SBarry Smith 
341a480d89SAdministrator    Notes: Click on the links below to see details on a particular solver
351a480d89SAdministrator 
361a480d89SAdministrator .seealso: PCSetType(), PC, PCCreate()
3776bdecfbSBarry Smith J*/
38a313700dSBarry Smith #define PCType char*
3982bf6240SBarry Smith #define PCNONE            "none"
4082bf6240SBarry Smith #define PCJACOBI          "jacobi"
4182bf6240SBarry Smith #define PCSOR             "sor"
4282bf6240SBarry Smith #define PCLU              "lu"
4382bf6240SBarry Smith #define PCSHELL           "shell"
4482bf6240SBarry Smith #define PCBJACOBI         "bjacobi"
4582bf6240SBarry Smith #define PCMG              "mg"
4682bf6240SBarry Smith #define PCEISENSTAT       "eisenstat"
4782bf6240SBarry Smith #define PCILU             "ilu"
4882bf6240SBarry Smith #define PCICC             "icc"
4982bf6240SBarry Smith #define PCASM             "asm"
50ab3e923fSDmitry Karpeev #define PCGASM            "gasm"
5194b7f48cSBarry Smith #define PCKSP             "ksp"
5282bf6240SBarry Smith #define PCCOMPOSITE       "composite"
53421c37bdSBarry Smith #define PCREDUNDANT       "redundant"
5427b520f0SBarry Smith #define PCSPAI            "spai"
55186905e3SBarry Smith #define PCNN              "nn"
564bbc92c1SBarry Smith #define PCCHOLESKY        "cholesky"
573a7fca6bSBarry Smith #define PCPBJACOBI        "pbjacobi"
587f5ff6fdSBarry Smith #define PCMAT             "mat"
59c4888f26SBarry Smith #define PCHYPRE           "hypre"
6037f80224SJose E. Roman #define PCPARMS           "parms"
610971522cSBarry Smith #define PCFIELDSPLIT      "fieldsplit"
62be16f70fSSatish Balay #define PCTFS             "tfs"
635582bec1SHong Zhang #define PCML              "ml"
6436a49750SBarry Smith #define PCPROMETHEUS      "prometheus"
652a6744ebSBarry Smith #define PCGALERKIN        "galerkin"
667233f9f0SBarry Smith #define PCEXOTIC          "exotic"
67bad7cb1dSBarry Smith #define PCHMPI            "hmpi"
68cf037197Ssdaitch #define PCSUPPORTGRAPH    "supportgraph"
69f4b8409dSBarry Smith #define PCASA             "asa"
7024c02a0fSBarry Smith #define PCCP              "cp"
71628b8437SMatthew Knepley #define PCBFBT            "bfbt"
72519d70e2SJed Brown #define PCLSC             "lsc"
731d6018f0SLisandro Dalcin #define PCPYTHON          "python"
74f91d8e95SBarry Smith #define PCPFMG            "pfmg"
75d851a50bSGlenn Hammond #define PCSYSPFMG         "syspfmg"
76df826632SBarry Smith #define PCREDISTRIBUTE    "redistribute"
77ac793be5SBarry Smith #define PCSVD             "svd"
78ac793be5SBarry Smith #define PCGAMG            "gamg"
79ac793be5SBarry Smith #define PCSACUSP          "sacusp"        /* these four run on NVIDIA GPUs using CUSP */
80bfc29b71SVictor Minden #define PCSACUSPPOLY      "sacusppoly"
818154be41SBarry Smith #define PCBICGSTABCUSP    "bicgstabcusp"
8204b59e88SVictor Minden #define PCAINVCUSP        "ainvcusp"
830c7d97c5SJed Brown #define PCBDDC            "bddc"
84123ea438SMatthew Knepley 
85123ea438SMatthew Knepley /* Logging support */
86*014dd563SJed Brown PETSC_EXTERN PetscClassId PC_CLASSID;
87123ea438SMatthew Knepley 
883d957683SBarry Smith /*E
893d957683SBarry Smith     PCSide - If the preconditioner is to be applied to the left, right
903d957683SBarry Smith      or symmetrically around the operator.
91d03aef70SBarry Smith 
923d957683SBarry Smith    Level: beginner
9321e95762SBarry Smith 
943d957683SBarry Smith .seealso:
953d957683SBarry Smith E*/
969e568555SJed Brown typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
979e568555SJed Brown #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
98*014dd563SJed Brown PETSC_EXTERN const char *PCSides[];
9972b7852fSLois Curfman McInnes 
100*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
101*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetType(PC,const PCType);
102*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetUp(PC);
103*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
104*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
105*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
106*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
107*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
108*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
109*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *);
110*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
1114d0a8057SBarry 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 
127*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
128*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *);
129*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool );
13084cb2905SBarry Smith 
131*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRegisterDestroy(void);
132*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRegisterAll(const char[]);
133*014dd563SJed Brown PETSC_EXTERN PetscBool PCRegisterAllCalled;
13484cb2905SBarry Smith 
135*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
13630de9b25SBarry Smith 
13730de9b25SBarry Smith /*MC
13830de9b25SBarry Smith    PCRegisterDynamic - Adds a method to the preconditioner package.
13930de9b25SBarry Smith 
14030de9b25SBarry Smith    Synopsis:
1411890ba74SBarry Smith    PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC))
14230de9b25SBarry Smith 
14330de9b25SBarry Smith    Not collective
14430de9b25SBarry Smith 
14530de9b25SBarry Smith    Input Parameters:
14630de9b25SBarry Smith +  name_solver - name of a new user-defined solver
14730de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
14830de9b25SBarry Smith .  name_create - name of routine to create method context
14930de9b25SBarry Smith -  routine_create - routine to create method context
15030de9b25SBarry Smith 
15130de9b25SBarry Smith    Notes:
15230de9b25SBarry Smith    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
15330de9b25SBarry Smith 
15430de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
15530de9b25SBarry Smith    is ignored.
15630de9b25SBarry Smith 
15730de9b25SBarry Smith    Sample usage:
15830de9b25SBarry Smith .vb
15930de9b25SBarry Smith    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
16030de9b25SBarry Smith               "MySolverCreate",MySolverCreate);
16130de9b25SBarry Smith .ve
16230de9b25SBarry Smith 
16330de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
16430de9b25SBarry Smith $     PCSetType(pc,"my_solver")
16530de9b25SBarry Smith    or at runtime via the option
16630de9b25SBarry Smith $     -pc_type my_solver
16730de9b25SBarry Smith 
16830de9b25SBarry Smith    Level: advanced
16930de9b25SBarry Smith 
170ab901514SBarry Smith    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
17130de9b25SBarry Smith            occuring in pathname will be replaced with appropriate values.
17230de9b25SBarry Smith          If your function is not being put into a shared library then use PCRegister() instead
17330de9b25SBarry Smith 
17430de9b25SBarry Smith .keywords: PC, register
17530de9b25SBarry Smith 
17630de9b25SBarry Smith .seealso: PCRegisterAll(), PCRegisterDestroy()
17730de9b25SBarry Smith M*/
178aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
179f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
1806df38c32SLois Curfman McInnes #else
181f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
1826df38c32SLois Curfman McInnes #endif
1836df38c32SLois Curfman McInnes 
184*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCReset(PC);
185*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
186*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);
187*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetType(PC,const PCType*);
18814c91fddSBarry Smith 
189*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
190*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
191*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
1925b116368SBarry Smith 
193*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
194*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
195*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
1964b0e389bSBarry Smith 
197*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer);
1987bc3d0afSSatish Balay 
199*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
200*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
201*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);
2028ed539a5SBarry Smith 
203*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
20471601f6fSBarry Smith 
205d6913704SBarry Smith /*
206d6913704SBarry Smith       These are used to provide extra scaling of preconditioned
2070f3b3ca1SHong Zhang    operator for time-stepping schemes like in SUNDIALS
208d6913704SBarry Smith */
209*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *);
210*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
211*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
212*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec);
213d6913704SBarry Smith 
21484cb2905SBarry Smith /* ------------- options specific to particular preconditioners --------- */
215329f5518SBarry Smith 
216*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
217*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC);
218*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC);
219*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
220*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
221*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
222d03aef70SBarry Smith 
223*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
224*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
225421c37bdSBarry Smith 
22615aa81f8SBarry Smith #define USE_PRECONDITIONER_MATRIX 0
22715aa81f8SBarry Smith #define USE_TRUE_MATRIX           1
228*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC);
229*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
230*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
2311eb62cbbSBarry Smith 
232*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPSetUseTrue(PC);
233981c4779SBarry Smith 
234*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
235*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
236*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
237*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
238*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
239*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
240*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
241*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**);
242*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*);
243*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
244*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]);
245aabeff55SBarry Smith 
246*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);
247d90ac83dSHong Zhang 
248*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType);
249*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal);
250d90ac83dSHong Zhang 
251*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
252*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
253*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC);
2542401956bSBarry Smith 
255*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
256*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal);
257*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
258421c37bdSBarry Smith 
259*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,const MatOrderingType);
260*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool );
261*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool );
262*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC);
263*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
264*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool );
265f5a88c2aSLois Curfman McInnes 
266*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
267*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
268b35a507dSBarry Smith 
269*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
270*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
271*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
272*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool );
273f746d493SDmitry Karpeev 
2743d957683SBarry Smith /*E
2753d957683SBarry Smith     PCASMType - Type of additive Schwarz method to use
2763d957683SBarry Smith 
277feb221f8SDmitry Karpeev $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
278feb221f8SDmitry Karpeev $                        and computed values in ghost regions are added together.
279feb221f8SDmitry Karpeev $                        Classical standard additive Schwarz.
280feb221f8SDmitry Karpeev $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
281feb221f8SDmitry Karpeev $                        region are discarded.
282feb221f8SDmitry Karpeev $                        Default.
283feb221f8SDmitry Karpeev $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
284feb221f8SDmitry Karpeev $                        region are added back in.
285feb221f8SDmitry Karpeev $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
286feb221f8SDmitry Karpeev $                        discarded.
287feb221f8SDmitry Karpeev $                        Not very good.
2883d957683SBarry Smith 
2893d957683SBarry Smith    Level: beginner
2903d957683SBarry Smith 
2913d957683SBarry Smith .seealso: PCASMSetType()
2923d957683SBarry Smith E*/
293d252947aSBarry Smith typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
294*014dd563SJed Brown PETSC_EXTERN const char *PCASMTypes[];
2959dcbbd2bSBarry Smith 
296*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
297*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
298*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
299*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
300*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
301*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
302981c4779SBarry Smith 
3033d957683SBarry Smith /*E
3046a4f0f73SDmitry Karpeev     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
305f746d493SDmitry Karpeev 
3066a4f0f73SDmitry Karpeev    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
3076a4f0f73SDmitry Karpeev    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
3086a4f0f73SDmitry Karpeev    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
3096a4f0f73SDmitry Karpeev    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
3106a4f0f73SDmitry Karpeev 
3116a4f0f73SDmitry Karpeev $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
3126a4f0f73SDmitry Karpeev $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
3136a4f0f73SDmitry Karpeev $                        from neighboring subdomains.
314feb221f8SDmitry Karpeev $                        Classical standard additive Schwarz.
3156a4f0f73SDmitry Karpeev $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
3166a4f0f73SDmitry Karpeev $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
3176a4f0f73SDmitry Karpeev $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
3186a4f0f73SDmitry Karpeev $                        assumption).
319feb221f8SDmitry Karpeev $                        Default.
3206a4f0f73SDmitry Karpeev $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
3216a4f0f73SDmitry Karpeev $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
3226a4f0f73SDmitry Karpeev $                        from neighboring subdomains.
3236a4f0f73SDmitry Karpeev $
3246a4f0f73SDmitry Karpeev $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
325feb221f8SDmitry Karpeev $                        Not very good.
326f746d493SDmitry Karpeev 
327f746d493SDmitry Karpeev    Level: beginner
328f746d493SDmitry Karpeev 
329f746d493SDmitry Karpeev .seealso: PCGASMSetType()
330f746d493SDmitry Karpeev E*/
331f746d493SDmitry Karpeev typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
332*014dd563SJed Brown PETSC_EXTERN const char *PCGASMTypes[];
333f746d493SDmitry Karpeev 
334*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
335*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool);
336*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
337*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );
338f746d493SDmitry Karpeev 
339*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
340*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]);
341*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
342*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
343*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
344*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);
345f746d493SDmitry Karpeev 
346f746d493SDmitry Karpeev /*E
3473d957683SBarry Smith     PCCompositeType - Determines how two or more preconditioner are composed
3483d957683SBarry Smith 
3493d957683SBarry Smith $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
3503d957683SBarry Smith $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
3513d957683SBarry Smith $                                computed after the previous preconditioner application
352ccb205f8SBarry Smith $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
353ccb205f8SBarry Smith $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
3543d957683SBarry Smith $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
3553d957683SBarry Smith $                         where first preconditioner is built from alpha I + S and second from
3563d957683SBarry Smith $                         alpha I + R
3573d957683SBarry Smith 
3583d957683SBarry Smith    Level: beginner
3593d957683SBarry Smith 
3603d957683SBarry Smith .seealso: PCCompositeSetType()
3613d957683SBarry Smith E*/
3623b224e63SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
363*014dd563SJed Brown PETSC_EXTERN const char *PCCompositeTypes[];
3649dcbbd2bSBarry Smith 
365*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSetUseTrue(PC);
366*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
367*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
368*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
369*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
370981c4779SBarry Smith 
371*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
372*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
373*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
374da3a660dSBarry Smith 
375*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
376*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
377*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
378*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
379*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
380*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
381*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
382*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
3833304466cSBarry Smith 
384*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
385*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
386*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
387*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
3883304466cSBarry Smith 
389*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
390*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
391*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
392*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
393*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
394*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
395084e4875SJed Brown 
396084e4875SJed Brown /*E
397084e4875SJed Brown     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
398084e4875SJed Brown 
399084e4875SJed Brown     Level: intermediate
400084e4875SJed Brown 
401084e4875SJed Brown .seealso: PCFieldSplitSchurPrecondition()
402084e4875SJed Brown E*/
403084e4875SJed Brown typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_DIAG,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
404*014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
405084e4875SJed Brown 
406ab1df9f5SJed Brown /*E
407c9c6ffaaSJed Brown     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
408ab1df9f5SJed Brown 
409ab1df9f5SJed Brown     Level: intermediate
410ab1df9f5SJed Brown 
411c9c6ffaaSJed Brown .seealso: PCFieldSplitSetSchurFactType()
412ab1df9f5SJed Brown E*/
413ab1df9f5SJed Brown typedef enum {
414c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
415c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
416c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
417c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_FULL
418c9c6ffaaSJed Brown } PCFieldSplitSchurFactType;
419*014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
420ab1df9f5SJed Brown 
421*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
422*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
423*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
4243d30b1ffSBarry Smith 
425*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
426*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
4272a6744ebSBarry Smith 
428*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
429*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);
4302102ba4dSSatish Balay 
431*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);
43267fac13cSBarry Smith 
433*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
434*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);
4356c699258SBarry Smith 
436*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
437*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);
4381b2093e4SBarry Smith 
439*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
440*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
441*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
4426c699258SBarry Smith 
443*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
444*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
445*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
446*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
44737f80224SJose E. Roman /*E
44837f80224SJose E. Roman     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
44937f80224SJose E. Roman 
45037f80224SJose E. Roman     Level: intermediate
45137f80224SJose E. Roman 
45237f80224SJose E. Roman .seealso: PCPARMSSetGlobal()
45337f80224SJose E. Roman E*/
45437f80224SJose E. Roman typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
455*014dd563SJed Brown PETSC_EXTERN const char *PCPARMSGlobalTypes[];
45637f80224SJose E. Roman /*E
45737f80224SJose E. Roman     PCPARMSLocalType - Determines the local preconditioner method in PARMS
45837f80224SJose E. Roman 
45937f80224SJose E. Roman     Level: intermediate
46037f80224SJose E. Roman 
46137f80224SJose E. Roman .seealso: PCPARMSSetLocal()
46237f80224SJose E. Roman E*/
46337f80224SJose E. Roman typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
464*014dd563SJed Brown PETSC_EXTERN const char *PCPARMSLocalTypes[];
46537f80224SJose E. Roman 
466*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type);
467*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type);
468*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits);
469*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart);
470*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym);
471*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2);
47237f80224SJose E. Roman 
473*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
474*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
475*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
476*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
477*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
478*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
479*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
48052284f78SMark F. Adams #define PCGAMGType char*
481*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,const PCGAMGType );
482*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n);
483*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n);
484*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool);
48556de70efSSatish Balay 
4860c7d97c5SJed Brown #if defined(PETSC_HAVE_PCBDDC)
48753cdbc3dSStefano Zampini /* Enum defining how to treat the coarse problem */
4880c7d97c5SJed Brown typedef enum {SEQUENTIAL_BDDC,REPLICATED_BDDC,PARALLEL_BDDC,MULTILEVEL_BDDC} CoarseProblemType;
489*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
490*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
491*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
492*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
493*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseProblemType(PC,CoarseProblemType);
494*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
4950c7d97c5SJed Brown #endif
4960c7d97c5SJed Brown 
497*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
498831a100dSStefano Zampini 
499d3b1e0e7SMatthew Knepley 
500123ea438SMatthew Knepley #endif /* __PETSCPC_H */
501