xref: /petsc/include/petscpc.h (revision efe48dd8cf8b935accbbb9f4bcb20bc83865fa4d)
1 /*
2       Preconditioner module.
3 */
4 #if !defined(__PETSCPC_H)
5 #define __PETSCPC_H
6 #include "petscmat.h"
7 #include "petscdm.h"
8 PETSC_EXTERN_CXX_BEGIN
9 
10 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT  PCInitializePackage(const char[]);
11 
12 /*
13     PCList contains the list of preconditioners currently registered
14    These are added with the PCRegisterDynamic() macro
15 */
16 extern PetscFList PCList;
17 
18 /*S
19      PC - Abstract PETSc object that manages all preconditioners
20 
21    Level: beginner
22 
23   Concepts: preconditioners
24 
25 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
26 S*/
27 typedef struct _p_PC* PC;
28 
29 /*E
30     PCType - String with the name of a PETSc preconditioner method or the creation function
31        with an optional dynamic library name, for example
32        http://www.mcs.anl.gov/petsc/lib.a:mypccreate()
33 
34    Level: beginner
35 
36    Notes: Click on the links below to see details on a particular solver
37 
38 .seealso: PCSetType(), PC, PCCreate()
39 E*/
40 #define PCType char*
41 #define PCNONE            "none"
42 #define PCJACOBI          "jacobi"
43 #define PCSOR             "sor"
44 #define PCLU              "lu"
45 #define PCSHELL           "shell"
46 #define PCBJACOBI         "bjacobi"
47 #define PCMG              "mg"
48 #define PCEISENSTAT       "eisenstat"
49 #define PCILU             "ilu"
50 #define PCICC             "icc"
51 #define PCASM             "asm"
52 #define PCGASM            "gasm"
53 #define PCKSP             "ksp"
54 #define PCCOMPOSITE       "composite"
55 #define PCREDUNDANT       "redundant"
56 #define PCSPAI            "spai"
57 #define PCNN              "nn"
58 #define PCCHOLESKY        "cholesky"
59 #define PCPBJACOBI        "pbjacobi"
60 #define PCMAT             "mat"
61 #define PCHYPRE           "hypre"
62 #define PCFIELDSPLIT      "fieldsplit"
63 #define PCTFS             "tfs"
64 #define PCML              "ml"
65 #define PCPROMETHEUS      "prometheus"
66 #define PCGALERKIN        "galerkin"
67 #define PCEXOTIC          "exotic"
68 #define PCOPENMP          "openmp"
69 #define PCSUPPORTGRAPH    "supportgraph"
70 #define PCASA             "asa"
71 #define PCCP              "cp"
72 #define PCBFBT            "bfbt"
73 #define PCLSC             "lsc"
74 #define PCPYTHON          "python"
75 #define PCPFMG            "pfmg"
76 #define PCSYSPFMG         "syspfmg"
77 #define PCREDISTRIBUTE    "redistribute"
78 #define PCSACUDA          "sacuda"
79 
80 /* Logging support */
81 extern PetscClassId PETSCKSP_DLLEXPORT PC_CLASSID;
82 
83 /*E
84     PCSide - If the preconditioner is to be applied to the left, right
85      or symmetrically around the operator.
86 
87    Level: beginner
88 
89 .seealso:
90 E*/
91 typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
92 extern const char *PCSides[];
93 
94 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCreate(MPI_Comm,PC*);
95 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetType(PC,const PCType);
96 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetUp(PC);
97 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetUpOnBlocks(PC);
98 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApply(PC,Vec,Vec);
99 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricLeft(PC,Vec,Vec);
100 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricRight(PC,Vec,Vec);
101 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
102 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyTranspose(PC,Vec,Vec);
103 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyTransposeExists(PC,PetscBool *);
104 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
105 
106 /*E
107     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
108 
109    Level: advanced
110 
111    Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
112 
113 .seealso: PCApplyRichardson()
114 E*/
115 typedef enum {
116               PCRICHARDSON_CONVERGED_RTOL               =  2,
117               PCRICHARDSON_CONVERGED_ATOL               =  3,
118               PCRICHARDSON_CONVERGED_ITS                =  4,
119               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
120 
121 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
122 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardsonExists(PC,PetscBool *);
123 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetInitialGuessNonzero(PC,PetscBool );
124 
125 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterDestroy(void);
126 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterAll(const char[]);
127 extern PetscBool  PCRegisterAllCalled;
128 
129 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
130 
131 /*MC
132    PCRegisterDynamic - Adds a method to the preconditioner package.
133 
134    Synopsis:
135    PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC))
136 
137    Not collective
138 
139    Input Parameters:
140 +  name_solver - name of a new user-defined solver
141 .  path - path (either absolute or relative) the library containing this solver
142 .  name_create - name of routine to create method context
143 -  routine_create - routine to create method context
144 
145    Notes:
146    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
147 
148    If dynamic libraries are used, then the fourth input argument (routine_create)
149    is ignored.
150 
151    Sample usage:
152 .vb
153    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
154               "MySolverCreate",MySolverCreate);
155 .ve
156 
157    Then, your solver can be chosen with the procedural interface via
158 $     PCSetType(pc,"my_solver")
159    or at runtime via the option
160 $     -pc_type my_solver
161 
162    Level: advanced
163 
164    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
165            occuring in pathname will be replaced with appropriate values.
166          If your function is not being put into a shared library then use PCRegister() instead
167 
168 .keywords: PC, register
169 
170 .seealso: PCRegisterAll(), PCRegisterDestroy()
171 M*/
172 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
173 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
174 #else
175 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
176 #endif
177 
178 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDestroy(PC);
179 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetFromOptions(PC);
180 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetType(PC,const PCType*);
181 
182 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatrix(PC,Mat*);
183 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
184 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
185 
186 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetOperators(PC,Mat,Mat,MatStructure);
187 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperators(PC,Mat*,Mat*,MatStructure*);
188 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
189 
190 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCView(PC,PetscViewer);
191 
192 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetOptionsPrefix(PC,const char[]);
193 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCAppendOptionsPrefix(PC,const char[]);
194 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOptionsPrefix(PC,const char*[]);
195 
196 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCComputeExplicitOperator(PC,Mat*);
197 
198 /*
199       These are used to provide extra scaling of preconditioned
200    operator for time-stepping schemes like in SUNDIALS
201 */
202 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetDiagonalScale(PC,PetscBool *);
203 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleLeft(PC,Vec,Vec);
204 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleRight(PC,Vec,Vec);
205 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetDiagonalScale(PC,Vec);
206 
207 /* ------------- options specific to particular preconditioners --------- */
208 
209 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowMax(PC);
210 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowSum(PC);
211 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseAbs(PC);
212 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric(PC,MatSORType);
213 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega(PC,PetscReal);
214 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations(PC,PetscInt,PetscInt);
215 
216 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatSetOmega(PC,PetscReal);
217 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatNoDiagonalScaling(PC);
218 
219 #define USE_PRECONDITIONER_MATRIX 0
220 #define USE_TRUE_MATRIX           1
221 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetUseTrueLocal(PC);
222 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
223 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
224 
225 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPSetUseTrue(PC);
226 
227 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
228 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
229 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
230 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
231 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
232 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
233 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
234 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetContext(PC,void**);
235 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetContext(PC,void*);
236 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetName(PC,const char[]);
237 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetName(PC,char*[]);
238 
239 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot(PC,PetscReal);
240 
241 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftType(PC,MatFactorShiftType);
242 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftAmount(PC,PetscReal);
243 
244 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
245 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
246 
247 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC,PetscReal);
248 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetColumnPivot(PC,PetscReal);
249 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
250 
251 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC,const MatOrderingType);
252 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC,PetscBool );
253 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC,PetscBool );
254 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC);
255 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC);
256 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC,PetscBool );
257 
258 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC,PetscInt);
259 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
260 
261 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
262 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
263 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetOverlap(PC,PetscInt);
264 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetSortIndices(PC,PetscBool );
265 
266 /*E
267     PCASMType - Type of additive Schwarz method to use
268 
269 $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
270 $                 and computed values in ghost regions are added together. Classical
271 $                 standard additive Schwarz
272 $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
273 $                    region are discarded. Default
274 $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
275 $                       region are added back in
276 $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
277 $                not very good.
278 
279    Level: beginner
280 
281 .seealso: PCASMSetType()
282 E*/
283 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
284 extern const char *PCASMTypes[];
285 
286 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetType(PC,PCASMType);
287 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
288 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMDestroySubdomains(PetscInt,IS[],IS[]);
289 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
290 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
291 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
292 
293 /*E
294     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per domain)
295 
296 $  PC_GASM_BASIC - symmetric version where residuals from the ghost points are used
297 $                 and computed values in ghost regions are added together. Classical
298 $                 standard additive Schwarz
299 $  PC_GASM_RESTRICT - residuals from ghost points are used but computed values in ghost
300 $                    region are discarded. Default
301 $  PC_GASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
302 $                       region are added back in
303 $  PC_GASM_NONE - ghost point residuals are not used, computed ghost values are discarded
304 $                not very good.
305 
306    Level: beginner
307 
308 .seealso: PCGASMSetType()
309 E*/
310 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
311 extern const char *PCGASMTypes[];
312 
313 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
314 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetTotalSubdomains(PC,PetscInt);
315 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetOverlap(PC,PetscInt);
316 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetSortIndices(PC,PetscBool );
317 
318 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetType(PC,PCGASMType);
319 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGASMCreateSubdomains(Mat,PetscInt,IS*[]);
320 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
321 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
322 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
323 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
324 
325 /*E
326     PCCompositeType - Determines how two or more preconditioner are composed
327 
328 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
329 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
330 $                                computed after the previous preconditioner application
331 $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
332 $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
333 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
334 $                         where first preconditioner is built from alpha I + S and second from
335 $                         alpha I + R
336 
337    Level: beginner
338 
339 .seealso: PCCompositeSetType()
340 E*/
341 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
342 extern const char *PCCompositeTypes[];
343 
344 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC);
345 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC,PCCompositeType);
346 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC,PCType);
347 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC,PetscInt,PC *);
348 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC,PetscScalar);
349 
350 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber(PC,PetscInt);
351 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter(PC,VecScatter,VecScatter);
352 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators(PC,Mat*,Mat*);
353 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC(PC,PC*);
354 
355 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetEpsilon(PC,double);
356 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetNBSteps(PC,PetscInt);
357 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMax(PC,PetscInt);
358 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMaxNew(PC,PetscInt);
359 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetBlockSize(PC,PetscInt);
360 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetCacheSize(PC,PetscInt);
361 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetVerbose(PC,PetscInt);
362 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetSp(PC,PetscInt);
363 
364 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC,const char[]);
365 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC,const char*[]);
366 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
367 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
368 
369 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*);
370 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC,PCCompositeType);
371 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC,PetscInt);
372 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC,const char[],IS);
373 
374 /*E
375     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
376 
377     Level: intermediate
378 
379 .seealso: PCFieldSplitSchurPrecondition()
380 E*/
381 typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_DIAG,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
382 extern const char *const PCFieldSplitSchurPreTypes[];
383 
384 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
385 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
386 
387 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetRestriction(PC,Mat);
388 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetInterpolation(PC,Mat);
389 
390 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetCoordinates(PC,PetscInt,PetscReal*);
391 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSASetVectors(PC,PetscInt,PetscReal *);
392 
393 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPythonSetType(PC,const char[]);
394 
395 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetDM(PC,DM);
396 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetDM(PC,DM*);
397 
398 
399 PETSC_EXTERN_CXX_END
400 
401 #endif /* __PETSCPC_H */
402