xref: /petsc/include/petscpc.h (revision cd47f5d92d6e538b4e65cbf3be7c5ae86cc86a49)
1d03aef70SBarry Smith /*
237f753daSBarry Smith       Preconditioner module.
3d03aef70SBarry Smith */
40a835dfdSSatish Balay #if !defined(__PETSCPC_H)
50a835dfdSSatish Balay #define __PETSCPC_H
60a835dfdSSatish Balay #include "petscmat.h"
7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
8d03aef70SBarry Smith 
9dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT  PCInitializePackage(const char[]);
101dbb0a54SBarry Smith 
11eec0b4cfSBarry Smith /*
12eec0b4cfSBarry Smith     PCList contains the list of preconditioners currently registered
13f1af5d2fSBarry Smith    These are added with the PCRegisterDynamic() macro
14eec0b4cfSBarry Smith */
15b0a32e0cSBarry Smith extern PetscFList PCList;
16f69a0ea3SMatthew Knepley #define PCType const char*
1782bf6240SBarry Smith 
183d957683SBarry Smith /*S
193d957683SBarry Smith      PC - Abstract PETSc object that manages all preconditioners
203d957683SBarry Smith 
213d957683SBarry Smith    Level: beginner
223d957683SBarry Smith 
233d957683SBarry Smith   Concepts: preconditioners
243d957683SBarry Smith 
251a480d89SAdministrator .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
263d957683SBarry Smith S*/
273d957683SBarry Smith typedef struct _p_PC* PC;
283d957683SBarry Smith 
293d957683SBarry Smith /*E
303d957683SBarry Smith     PCType - String with the name of a PETSc preconditioner method or the creation function
313d957683SBarry Smith        with an optional dynamic library name, for example
323d957683SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mypccreate()
333d957683SBarry Smith 
343d957683SBarry Smith    Level: beginner
353d957683SBarry Smith 
361a480d89SAdministrator    Notes: Click on the links below to see details on a particular solver
371a480d89SAdministrator 
381a480d89SAdministrator .seealso: PCSetType(), PC, PCCreate()
393d957683SBarry Smith E*/
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"
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"
57a6b93e97SBarry Smith #define PCSAMG            "samg"
583a7fca6bSBarry Smith #define PCPBJACOBI        "pbjacobi"
597f5ff6fdSBarry Smith #define PCMAT             "mat"
60c4888f26SBarry Smith #define PCHYPRE           "hypre"
610971522cSBarry Smith #define PCFIELDSPLIT      "fieldsplit"
62be16f70fSSatish Balay #define PCTFS             "tfs"
635582bec1SHong Zhang #define PCML              "ml"
6436a49750SBarry Smith #define PCPROMETHEUS      "prometheus"
65123ea438SMatthew Knepley 
66123ea438SMatthew Knepley /* Logging support */
67dba47a55SKris Buschelman extern PetscCookie PETSCKSP_DLLEXPORT PC_COOKIE;
68123ea438SMatthew Knepley 
693d957683SBarry Smith /*E
703d957683SBarry Smith     PCSide - If the preconditioner is to be applied to the left, right
713d957683SBarry Smith      or symmetrically around the operator.
72d03aef70SBarry Smith 
733d957683SBarry Smith    Level: beginner
7421e95762SBarry Smith 
753d957683SBarry Smith .seealso:
763d957683SBarry Smith E*/
77aad2872bSLois Curfman McInnes typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
789dcbbd2bSBarry Smith extern const char *PCSides[];
7972b7852fSLois Curfman McInnes 
80dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCreate(MPI_Comm,PC*);
81f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetType(PC,PCType);
82dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetUp(PC);
83dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetUpOnBlocks(PC);
84dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApply(PC,Vec,Vec);
85dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricLeft(PC,Vec,Vec);
86dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricRight(PC,Vec,Vec);
87dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
88dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyTranspose(PC,Vec,Vec);
89dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCHasApplyTranspose(PC,PetscTruth*);
90dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
91dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt);
92dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardsonExists(PC,PetscTruth*);
9384cb2905SBarry Smith 
94dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterDestroy(void);
95dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterAll(const char[]);
962bad1931SBarry Smith extern PetscTruth PCRegisterAllCalled;
9784cb2905SBarry Smith 
98dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
9930de9b25SBarry Smith 
10030de9b25SBarry Smith /*MC
10130de9b25SBarry Smith    PCRegisterDynamic - Adds a method to the preconditioner package.
10230de9b25SBarry Smith 
10330de9b25SBarry Smith    Synopsis:
104d360dc6fSBarry Smith    PetscErrorCode PCRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(PC))
10530de9b25SBarry Smith 
10630de9b25SBarry Smith    Not collective
10730de9b25SBarry Smith 
10830de9b25SBarry Smith    Input Parameters:
10930de9b25SBarry Smith +  name_solver - name of a new user-defined solver
11030de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
11130de9b25SBarry Smith .  name_create - name of routine to create method context
11230de9b25SBarry Smith -  routine_create - routine to create method context
11330de9b25SBarry Smith 
11430de9b25SBarry Smith    Notes:
11530de9b25SBarry Smith    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
11630de9b25SBarry Smith 
11730de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
11830de9b25SBarry Smith    is ignored.
11930de9b25SBarry Smith 
12030de9b25SBarry Smith    Sample usage:
12130de9b25SBarry Smith .vb
12230de9b25SBarry Smith    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
12330de9b25SBarry Smith               "MySolverCreate",MySolverCreate);
12430de9b25SBarry Smith .ve
12530de9b25SBarry Smith 
12630de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
12730de9b25SBarry Smith $     PCSetType(pc,"my_solver")
12830de9b25SBarry Smith    or at runtime via the option
12930de9b25SBarry Smith $     -pc_type my_solver
13030de9b25SBarry Smith 
13130de9b25SBarry Smith    Level: advanced
13230de9b25SBarry Smith 
133ab901514SBarry Smith    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
13430de9b25SBarry Smith            occuring in pathname will be replaced with appropriate values.
13530de9b25SBarry Smith          If your function is not being put into a shared library then use PCRegister() instead
13630de9b25SBarry Smith 
13730de9b25SBarry Smith .keywords: PC, register
13830de9b25SBarry Smith 
13930de9b25SBarry Smith .seealso: PCRegisterAll(), PCRegisterDestroy()
14030de9b25SBarry Smith M*/
141aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
142f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
1436df38c32SLois Curfman McInnes #else
144f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
1456df38c32SLois Curfman McInnes #endif
1466df38c32SLois Curfman McInnes 
147dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDestroy(PC);
148dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetFromOptions(PC);
149dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetType(PC,PCType*);
15014c91fddSBarry Smith 
151dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetFactoredMatrix(PC,Mat*);
152dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
153dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
1545b116368SBarry Smith 
155dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetOperators(PC,Mat,Mat,MatStructure);
156dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperators(PC,Mat*,Mat*,MatStructure*);
1574b0e389bSBarry Smith 
158dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCView(PC,PetscViewer);
1597bc3d0afSSatish Balay 
160dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetOptionsPrefix(PC,const char[]);
161dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCAppendOptionsPrefix(PC,const char[]);
1622dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOptionsPrefix(PC,const char*[]);
1638ed539a5SBarry Smith 
164dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCComputeExplicitOperator(PC,Mat*);
16571601f6fSBarry Smith 
166d6913704SBarry Smith /*
167d6913704SBarry Smith       These are used to provide extra scaling of preconditioned
1680f3b3ca1SHong Zhang    operator for time-stepping schemes like in SUNDIALS
169d6913704SBarry Smith */
170dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScale(PC,PetscTruth*);
171dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleLeft(PC,Vec,Vec);
172dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleRight(PC,Vec,Vec);
173dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleSet(PC,Vec);
174d6913704SBarry Smith 
17584cb2905SBarry Smith /* ------------- options specific to particular preconditioners --------- */
176329f5518SBarry Smith 
177dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowMax(PC);
178*cd47f5d9SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseAbs(PC);
179dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric(PC,MatSORType);
180dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega(PC,PetscReal);
181dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations(PC,PetscInt,PetscInt);
182d03aef70SBarry Smith 
183dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatSetOmega(PC,PetscReal);
184dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatNoDiagonalScaling(PC);
185421c37bdSBarry Smith 
18615aa81f8SBarry Smith #define USE_PRECONDITIONER_MATRIX 0
18715aa81f8SBarry Smith #define USE_TRUE_MATRIX           1
188dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetUseTrueLocal(PC);
189dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
190dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
1911eb62cbbSBarry Smith 
192dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPSetUseTrue(PC);
193981c4779SBarry Smith 
194be29d3c6SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApply(PC,PetscErrorCode (*)(void*,Vec,Vec));
195dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyTranspose(PC,PetscErrorCode (*)(void*,Vec,Vec));
196dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetSetUp(PC,PetscErrorCode (*)(void*));
197be29d3c6SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyRichardson(PC,PetscErrorCode (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt));
198dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetView(PC,PetscErrorCode (*)(void*,PetscViewer));
19918be62a5SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetDestroy(PC,PetscErrorCode (*)(void*));
200be29d3c6SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetContext(PC,void**);
201be29d3c6SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetContext(PC,void*);
202dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetName(PC,const char[]);
203dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetName(PC,char*[]);
204aabeff55SBarry Smith 
205dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot(PC,PetscReal);
206dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero(PC,PetscReal);
207dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd(PC,PetscTruth);
208ee45ca4aSHong Zhang 
2092401956bSBarry Smith 
21055ba2a51SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC,PetscReal);
2112401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC,PetscReal);
2122401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
213421c37bdSBarry Smith 
2142401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrdering(PC,MatOrderingType);
2152401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC,PetscTruth);
2162401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC,PetscTruth);
2172401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC);
21812a54f56SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC);
2192401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC,PetscTruth);
220f5a88c2aSLois Curfman McInnes 
2212401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC,PetscInt);
2222401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC,PetscReal,PetscReal,PetscInt);
223b35a507dSBarry Smith 
224dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetLocalSubdomains(PC,PetscInt,IS[]);
225dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetTotalSubdomains(PC,PetscInt,IS[]);
226dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetOverlap(PC,PetscInt);
2273d957683SBarry Smith /*E
2283d957683SBarry Smith     PCASMType - Type of additive Schwarz method to use
2293d957683SBarry Smith 
2303d957683SBarry Smith $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
2313d957683SBarry Smith $                 and computed values in ghost regions are added together. Classical
2323d957683SBarry Smith $                 standard additive Schwarz
2333d957683SBarry Smith $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
2343d957683SBarry Smith $                    region are discarded. Default
2353d957683SBarry Smith $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
2363d957683SBarry Smith $                       region are added back in
2373d957683SBarry Smith $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
2383d957683SBarry Smith $                not very good.
2393d957683SBarry Smith 
2403d957683SBarry Smith    Level: beginner
2413d957683SBarry Smith 
2423d957683SBarry Smith .seealso: PCASMSetType()
2433d957683SBarry Smith E*/
244d252947aSBarry Smith typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
2459dcbbd2bSBarry Smith extern const char *PCASMTypes[];
2469dcbbd2bSBarry Smith 
247dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetType(PC,PCASMType);
248dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt *,IS **);
249dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetUseInPlace(PC);
250dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubdomains(PC,PetscInt*,IS*[]);
251dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
252981c4779SBarry Smith 
2533d957683SBarry Smith /*E
2543d957683SBarry Smith     PCCompositeType - Determines how two or more preconditioner are composed
2553d957683SBarry Smith 
2563d957683SBarry Smith $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
2573d957683SBarry Smith $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
2583d957683SBarry Smith $                                computed after the previous preconditioner application
2593d957683SBarry Smith $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
2603d957683SBarry Smith $                         where first preconditioner is built from alpha I + S and second from
2613d957683SBarry Smith $                         alpha I + R
2623d957683SBarry Smith 
2633d957683SBarry Smith    Level: beginner
2643d957683SBarry Smith 
2653d957683SBarry Smith .seealso: PCCompositeSetType()
2663d957683SBarry Smith E*/
2671d1367b7SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType;
2689dcbbd2bSBarry Smith extern const char *PCCompositeTypes[];
2699dcbbd2bSBarry Smith 
270dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC);
271dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC,PCCompositeType);
272dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC,PCType);
273dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC pc,PetscInt n,PC *);
274dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC,PetscScalar);
275981c4779SBarry Smith 
276dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter(PC,VecScatter,VecScatter);
277dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators(PC,Mat*,Mat*);
278dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC(PC,PC*);
279da3a660dSBarry Smith 
280dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetEpsilon(PC,double);
281dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetNBSteps(PC,PetscInt);
282dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMax(PC,PetscInt);
283dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMaxNew(PC,PetscInt);
284dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetBlockSize(PC,PetscInt);
285dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetCacheSize(PC,PetscInt);
286dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetVerbose(PC,PetscInt);
287dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetSp(PC,PetscInt);
2883304466cSBarry Smith 
289dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC,const char[]);
290dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
291dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
2923304466cSBarry Smith 
293dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC,PetscInt,PetscInt*);
294dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC,PCCompositeType);
2953d30b1ffSBarry Smith 
296c02b24c5SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetCoordinates(PC,PetscReal*);
29767fac13cSBarry Smith 
298e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
299123ea438SMatthew Knepley #endif /* __PETSCPC_H */
300