xref: /petsc/include/petscpc.h (revision d7d82daadd7e0d6eee5c84c3bd60dd428ccf339b)
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 
9dfbe8321SBarry Smith EXTERN PetscErrorCode 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;
1649773a63SBarry Smith #define PCType 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"
57237f4ff4SBarry Smith #define PCRAMG       "ramg"
58a6b93e97SBarry Smith #define PCSAMG       "samg"
593a7fca6bSBarry Smith #define PCPBJACOBI   "pbjacobi"
607f5ff6fdSBarry Smith #define PCMAT        "mat"
61c4888f26SBarry Smith #define PCHYPRE      "hypre"
62123ea438SMatthew Knepley 
63123ea438SMatthew Knepley /* Logging support */
646849ba73SBarry Smith extern PetscCookie PC_COOKIE;
656849ba73SBarry Smith extern PetscEvent    PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
666849ba73SBarry Smith extern PetscEvent    PC_ApplySymmetricRight, PC_ModifySubMatrices;
67123ea438SMatthew Knepley 
683d957683SBarry Smith /*E
693d957683SBarry Smith     PCSide - If the preconditioner is to be applied to the left, right
703d957683SBarry Smith      or symmetrically around the operator.
71d03aef70SBarry Smith 
723d957683SBarry Smith    Level: beginner
7321e95762SBarry Smith 
743d957683SBarry Smith .seealso:
753d957683SBarry Smith E*/
76aad2872bSLois Curfman McInnes typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
7772b7852fSLois Curfman McInnes 
78dfbe8321SBarry Smith EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
79dfbe8321SBarry Smith EXTERN PetscErrorCode PCSetType(PC,const PCType);
80dfbe8321SBarry Smith EXTERN PetscErrorCode PCSetUp(PC);
81dfbe8321SBarry Smith EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
82dfbe8321SBarry Smith EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
83dfbe8321SBarry Smith EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
84dfbe8321SBarry Smith EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
85dfbe8321SBarry Smith EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
86dfbe8321SBarry Smith EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
87dfbe8321SBarry Smith EXTERN PetscErrorCode PCHasApplyTranspose(PC,PetscTruth*);
88dfbe8321SBarry Smith EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
89dfbe8321SBarry Smith EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,int);
90dfbe8321SBarry Smith EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscTruth*);
9184cb2905SBarry Smith 
92dfbe8321SBarry Smith EXTERN PetscErrorCode        PCRegisterDestroy(void);
93dfbe8321SBarry Smith EXTERN PetscErrorCode        PCRegisterAll(const char[]);
942bad1931SBarry Smith extern PetscTruth PCRegisterAllCalled;
9584cb2905SBarry Smith 
966849ba73SBarry Smith EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
9730de9b25SBarry Smith 
9830de9b25SBarry Smith /*MC
9930de9b25SBarry Smith    PCRegisterDynamic - Adds a method to the preconditioner package.
10030de9b25SBarry Smith 
10130de9b25SBarry Smith    Synopsis:
1026849ba73SBarry Smith    int PCRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(PC))
10330de9b25SBarry Smith 
10430de9b25SBarry Smith    Not collective
10530de9b25SBarry Smith 
10630de9b25SBarry Smith    Input Parameters:
10730de9b25SBarry Smith +  name_solver - name of a new user-defined solver
10830de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
10930de9b25SBarry Smith .  name_create - name of routine to create method context
11030de9b25SBarry Smith -  routine_create - routine to create method context
11130de9b25SBarry Smith 
11230de9b25SBarry Smith    Notes:
11330de9b25SBarry Smith    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
11430de9b25SBarry Smith 
11530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
11630de9b25SBarry Smith    is ignored.
11730de9b25SBarry Smith 
11830de9b25SBarry Smith    Sample usage:
11930de9b25SBarry Smith .vb
12030de9b25SBarry Smith    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
12130de9b25SBarry Smith               "MySolverCreate",MySolverCreate);
12230de9b25SBarry Smith .ve
12330de9b25SBarry Smith 
12430de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
12530de9b25SBarry Smith $     PCSetType(pc,"my_solver")
12630de9b25SBarry Smith    or at runtime via the option
12730de9b25SBarry Smith $     -pc_type my_solver
12830de9b25SBarry Smith 
12930de9b25SBarry Smith    Level: advanced
13030de9b25SBarry Smith 
13130de9b25SBarry Smith    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT}, or ${any environmental variable}
13230de9b25SBarry Smith            occuring in pathname will be replaced with appropriate values.
13330de9b25SBarry Smith          If your function is not being put into a shared library then use PCRegister() instead
13430de9b25SBarry Smith 
13530de9b25SBarry Smith .keywords: PC, register
13630de9b25SBarry Smith 
13730de9b25SBarry Smith .seealso: PCRegisterAll(), PCRegisterDestroy()
13830de9b25SBarry Smith M*/
139aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
140f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
1416df38c32SLois Curfman McInnes #else
142f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
1436df38c32SLois Curfman McInnes #endif
1446df38c32SLois Curfman McInnes 
145dfbe8321SBarry Smith EXTERN PetscErrorCode PCDestroy(PC);
146dfbe8321SBarry Smith EXTERN PetscErrorCode PCSetFromOptions(PC);
147dfbe8321SBarry Smith EXTERN PetscErrorCode PCGetType(PC,PCType*);
14814c91fddSBarry Smith 
149dfbe8321SBarry Smith EXTERN PetscErrorCode PCGetFactoredMatrix(PC,Mat*);
150dfbe8321SBarry Smith EXTERN PetscErrorCode PCSetModifySubMatrices(PC,int(*)(PC,int,const IS[],const IS[],Mat[],void*),void*);
151dfbe8321SBarry Smith EXTERN PetscErrorCode PCModifySubMatrices(PC,int,const IS[],const IS[],Mat[],void*);
1525b116368SBarry Smith 
153dfbe8321SBarry Smith EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
154dfbe8321SBarry Smith EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
1554b0e389bSBarry Smith 
156dfbe8321SBarry Smith EXTERN PetscErrorCode PCView(PC,PetscViewer);
1577bc3d0afSSatish Balay 
158dfbe8321SBarry Smith EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
159dfbe8321SBarry Smith EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
160dfbe8321SBarry Smith EXTERN PetscErrorCode PCGetOptionsPrefix(PC,char*[]);
1618ed539a5SBarry Smith 
162dfbe8321SBarry Smith EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
16371601f6fSBarry Smith 
164d6913704SBarry Smith /*
165d6913704SBarry Smith       These are used to provide extra scaling of preconditioned
166d6913704SBarry Smith    operator for time-stepping schemes like in PVODE
167d6913704SBarry Smith */
168dfbe8321SBarry Smith EXTERN PetscErrorCode PCDiagonalScale(PC,PetscTruth*);
169dfbe8321SBarry Smith EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
170dfbe8321SBarry Smith EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
171dfbe8321SBarry Smith EXTERN PetscErrorCode PCDiagonalScaleSet(PC,Vec);
172d6913704SBarry Smith 
17384cb2905SBarry Smith /* ------------- options specific to particular preconditioners --------- */
174329f5518SBarry Smith 
175dfbe8321SBarry Smith EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
176dfbe8321SBarry Smith EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
177dfbe8321SBarry Smith EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
178dfbe8321SBarry Smith EXTERN PetscErrorCode PCSORSetIterations(PC,int,int);
179d03aef70SBarry Smith 
180dfbe8321SBarry Smith EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
181dfbe8321SBarry Smith EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
182421c37bdSBarry Smith 
18315aa81f8SBarry Smith #define USE_PRECONDITIONER_MATRIX 0
18415aa81f8SBarry Smith #define USE_TRUE_MATRIX           1
185dfbe8321SBarry Smith EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC);
186dfbe8321SBarry Smith EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,int,const int[]);
187dfbe8321SBarry Smith EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,int,const int[]);
1881eb62cbbSBarry Smith 
189dfbe8321SBarry Smith EXTERN PetscErrorCode PCKSPSetUseTrue(PC);
190981c4779SBarry Smith 
1916849ba73SBarry Smith EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(void*,Vec,Vec),void*);
1926849ba73SBarry Smith EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(void*,Vec,Vec));
1936849ba73SBarry Smith EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(void*));
1946849ba73SBarry Smith EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,int),void*);
1956849ba73SBarry Smith EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(void*,PetscViewer));
196dfbe8321SBarry Smith EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
197dfbe8321SBarry Smith EXTERN PetscErrorCode PCShellGetName(PC,char*[]);
198aabeff55SBarry Smith 
199dfbe8321SBarry Smith EXTERN PetscErrorCode PCLUSetMatOrdering(PC,MatOrderingType);
200dfbe8321SBarry Smith EXTERN PetscErrorCode PCLUSetReuseOrdering(PC,PetscTruth);
201dfbe8321SBarry Smith EXTERN PetscErrorCode PCLUSetReuseFill(PC,PetscTruth);
202dfbe8321SBarry Smith EXTERN PetscErrorCode PCLUSetUseInPlace(PC);
203dfbe8321SBarry Smith EXTERN PetscErrorCode PCLUSetFill(PC,PetscReal);
204dfbe8321SBarry Smith EXTERN PetscErrorCode PCLUSetDamping(PC,PetscReal);
205dfbe8321SBarry Smith EXTERN PetscErrorCode PCLUSetShift(PC,PetscTruth);
206dfbe8321SBarry Smith EXTERN PetscErrorCode PCLUSetPivoting(PC,PetscReal);
207dfbe8321SBarry Smith EXTERN PetscErrorCode PCLUSetPivotInBlocks(PC,PetscTruth);
208dfbe8321SBarry Smith EXTERN PetscErrorCode PCLUSetZeroPivot(PC,PetscReal);
209421c37bdSBarry Smith 
210dfbe8321SBarry Smith EXTERN PetscErrorCode PCCholeskySetMatOrdering(PC,MatOrderingType);
211dfbe8321SBarry Smith EXTERN PetscErrorCode PCCholeskySetReuseOrdering(PC,PetscTruth);
212dfbe8321SBarry Smith EXTERN PetscErrorCode PCCholeskySetReuseFill(PC,PetscTruth);
213dfbe8321SBarry Smith EXTERN PetscErrorCode PCCholeskySetUseInPlace(PC);
214dfbe8321SBarry Smith EXTERN PetscErrorCode PCCholeskySetFill(PC,PetscReal);
215dfbe8321SBarry Smith EXTERN PetscErrorCode PCCholeskySetDamping(PC,PetscReal);
216dfbe8321SBarry Smith EXTERN PetscErrorCode PCCholeskySetShift(PC,PetscTruth);
217dfbe8321SBarry Smith EXTERN PetscErrorCode PCCholeskySetPivotInBlocks(PC,PetscTruth);
218f5a88c2aSLois Curfman McInnes 
219dfbe8321SBarry Smith EXTERN PetscErrorCode PCILUSetMatOrdering(PC,MatOrderingType);
220dfbe8321SBarry Smith EXTERN PetscErrorCode PCILUSetUseInPlace(PC);
221dfbe8321SBarry Smith EXTERN PetscErrorCode PCILUSetFill(PC,PetscReal);
222*d7d82daaSBarry Smith EXTERN PetscErrorCode PCILUSetLevels(PC,PetscInt);
223dfbe8321SBarry Smith EXTERN PetscErrorCode PCILUSetReuseOrdering(PC,PetscTruth);
224*d7d82daaSBarry Smith EXTERN PetscErrorCode PCILUSetUseDropTolerance(PC,PetscReal,PetscReal,PetscInt);
225dfbe8321SBarry Smith EXTERN PetscErrorCode PCILUDTSetReuseFill(PC,PetscTruth);
226dfbe8321SBarry Smith EXTERN PetscErrorCode PCILUSetAllowDiagonalFill(PC);
227dfbe8321SBarry Smith EXTERN PetscErrorCode PCILUSetDamping(PC,PetscReal);
228dfbe8321SBarry Smith EXTERN PetscErrorCode PCILUSetShift(PC,PetscTruth);
229dfbe8321SBarry Smith EXTERN PetscErrorCode PCILUSetPivotInBlocks(PC,PetscTruth);
230dfbe8321SBarry Smith EXTERN PetscErrorCode PCILUSetZeroPivot(PC,PetscReal);
231a2ce50c7SBarry Smith 
232dfbe8321SBarry Smith EXTERN PetscErrorCode PCICCSetMatOrdering(PC,MatOrderingType);
233dfbe8321SBarry Smith EXTERN PetscErrorCode PCICCSetFill(PC,PetscReal);
234*d7d82daaSBarry Smith EXTERN PetscErrorCode PCICCSetLevels(PC,PetscInt);
235dfbe8321SBarry Smith EXTERN PetscErrorCode PCICCSetDamping(PC,PetscReal);
236dfbe8321SBarry Smith EXTERN PetscErrorCode PCICCSetShift(PC,PetscTruth);
237dfbe8321SBarry Smith EXTERN PetscErrorCode PCICCSetPivotInBlocks(PC,PetscTruth);
238dfbe8321SBarry Smith EXTERN PetscErrorCode PCICCSetZeroPivot(PC,PetscReal);
239b35a507dSBarry Smith 
240*d7d82daaSBarry Smith EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[]);
241*d7d82daaSBarry Smith EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[]);
242*d7d82daaSBarry Smith EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
2433d957683SBarry Smith /*E
2443d957683SBarry Smith     PCASMType - Type of additive Schwarz method to use
2453d957683SBarry Smith 
2463d957683SBarry Smith $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
2473d957683SBarry Smith $                 and computed values in ghost regions are added together. Classical
2483d957683SBarry Smith $                 standard additive Schwarz
2493d957683SBarry Smith $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
2503d957683SBarry Smith $                    region are discarded. Default
2513d957683SBarry Smith $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
2523d957683SBarry Smith $                       region are added back in
2533d957683SBarry Smith $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
2543d957683SBarry Smith $                not very good.
2553d957683SBarry Smith 
2563d957683SBarry Smith    Level: beginner
2573d957683SBarry Smith 
2583d957683SBarry Smith .seealso: PCASMSetType()
2593d957683SBarry Smith E*/
260d252947aSBarry Smith typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
261dfbe8321SBarry Smith EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
262dfbe8321SBarry Smith EXTERN PetscErrorCode PCASMCreateSubdomains2D(int,int,int,int,int,int,int *,IS **);
263dfbe8321SBarry Smith EXTERN PetscErrorCode PCASMSetUseInPlace(PC);
264dfbe8321SBarry Smith EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,int*,IS*[]);
265dfbe8321SBarry Smith EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,int*,Mat*[]);
266981c4779SBarry Smith 
2673d957683SBarry Smith /*E
2683d957683SBarry Smith     PCCompositeType - Determines how two or more preconditioner are composed
2693d957683SBarry Smith 
2703d957683SBarry Smith $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
2713d957683SBarry Smith $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
2723d957683SBarry Smith $                                computed after the previous preconditioner application
2733d957683SBarry Smith $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
2743d957683SBarry Smith $                         where first preconditioner is built from alpha I + S and second from
2753d957683SBarry Smith $                         alpha I + R
2763d957683SBarry Smith 
2773d957683SBarry Smith    Level: beginner
2783d957683SBarry Smith 
2793d957683SBarry Smith .seealso: PCCompositeSetType()
2803d957683SBarry Smith E*/
2811d1367b7SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType;
282dfbe8321SBarry Smith EXTERN PetscErrorCode PCCompositeSetUseTrue(PC);
283dfbe8321SBarry Smith EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
284dfbe8321SBarry Smith EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
285dfbe8321SBarry Smith EXTERN PetscErrorCode PCCompositeGetPC(PC pc,int n,PC *);
286dfbe8321SBarry Smith EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
287981c4779SBarry Smith 
288dfbe8321SBarry Smith EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
289dfbe8321SBarry Smith EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
290dfbe8321SBarry Smith EXTERN PetscErrorCode PCRedundantGetPC(PC,PC*);
291dfbe8321SBarry Smith EXTERN PetscErrorCode MatGetOrderingList(PetscFList *list);
292da3a660dSBarry Smith 
293dfbe8321SBarry Smith EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
294dfbe8321SBarry Smith EXTERN PetscErrorCode PCSPAISetNBSteps(PC,int);
295dfbe8321SBarry Smith EXTERN PetscErrorCode PCSPAISetMax(PC,int);
296dfbe8321SBarry Smith EXTERN PetscErrorCode PCSPAISetMaxNew(PC,int);
297dfbe8321SBarry Smith EXTERN PetscErrorCode PCSPAISetBlockSize(PC,int);
298dfbe8321SBarry Smith EXTERN PetscErrorCode PCSPAISetCacheSize(PC,int);
299dfbe8321SBarry Smith EXTERN PetscErrorCode PCSPAISetVerbose(PC,int);
300dfbe8321SBarry Smith EXTERN PetscErrorCode PCSPAISetSp(PC,int);
3013304466cSBarry Smith 
302dfbe8321SBarry Smith EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
303dfbe8321SBarry Smith EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,int*,const int*[]);
304dfbe8321SBarry Smith EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,int*,const int*[]);
3053304466cSBarry Smith 
306e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
307123ea438SMatthew Knepley #endif /* __PETSCPC_H */
308