xref: /petsc/include/petscpc.h (revision 09573ac72a50d3e7ecd55a2b7f0ef28450cd0a8b)
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 #define PCSVD             "svd"
80 
81 /* Logging support */
82 extern PetscClassId PETSCKSP_DLLEXPORT PC_CLASSID;
83 
84 /*E
85     PCSide - If the preconditioner is to be applied to the left, right
86      or symmetrically around the operator.
87 
88    Level: beginner
89 
90 .seealso:
91 E*/
92 typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
93 extern const char *PCSides[];
94 
95 extern PetscErrorCode PETSCKSP_DLLEXPORT PCCreate(MPI_Comm,PC*);
96 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSetType(PC,const PCType);
97 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSetUp(PC);
98 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSetUpOnBlocks(PC);
99 extern PetscErrorCode PETSCKSP_DLLEXPORT PCApply(PC,Vec,Vec);
100 extern PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricLeft(PC,Vec,Vec);
101 extern PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricRight(PC,Vec,Vec);
102 extern PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
103 extern PetscErrorCode PETSCKSP_DLLEXPORT PCApplyTranspose(PC,Vec,Vec);
104 extern PetscErrorCode PETSCKSP_DLLEXPORT PCApplyTransposeExists(PC,PetscBool *);
105 extern PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
106 
107 /*E
108     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
109 
110    Level: advanced
111 
112    Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
113 
114 .seealso: PCApplyRichardson()
115 E*/
116 typedef enum {
117               PCRICHARDSON_CONVERGED_RTOL               =  2,
118               PCRICHARDSON_CONVERGED_ATOL               =  3,
119               PCRICHARDSON_CONVERGED_ITS                =  4,
120               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
121 
122 extern PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
123 extern PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardsonExists(PC,PetscBool *);
124 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSetInitialGuessNonzero(PC,PetscBool );
125 
126 extern PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterDestroy(void);
127 extern PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterAll(const char[]);
128 extern PetscBool  PCRegisterAllCalled;
129 
130 extern PetscErrorCode PETSCKSP_DLLEXPORT PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
131 
132 /*MC
133    PCRegisterDynamic - Adds a method to the preconditioner package.
134 
135    Synopsis:
136    PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC))
137 
138    Not collective
139 
140    Input Parameters:
141 +  name_solver - name of a new user-defined solver
142 .  path - path (either absolute or relative) the library containing this solver
143 .  name_create - name of routine to create method context
144 -  routine_create - routine to create method context
145 
146    Notes:
147    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
148 
149    If dynamic libraries are used, then the fourth input argument (routine_create)
150    is ignored.
151 
152    Sample usage:
153 .vb
154    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
155               "MySolverCreate",MySolverCreate);
156 .ve
157 
158    Then, your solver can be chosen with the procedural interface via
159 $     PCSetType(pc,"my_solver")
160    or at runtime via the option
161 $     -pc_type my_solver
162 
163    Level: advanced
164 
165    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
166            occuring in pathname will be replaced with appropriate values.
167          If your function is not being put into a shared library then use PCRegister() instead
168 
169 .keywords: PC, register
170 
171 .seealso: PCRegisterAll(), PCRegisterDestroy()
172 M*/
173 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
174 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
175 #else
176 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
177 #endif
178 
179 extern PetscErrorCode PETSCKSP_DLLEXPORT PCDestroy(PC);
180 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSetFromOptions(PC);
181 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGetType(PC,const PCType*);
182 
183 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatrix(PC,Mat*);
184 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
185 extern PetscErrorCode PETSCKSP_DLLEXPORT PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
186 
187 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSetOperators(PC,Mat,Mat,MatStructure);
188 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperators(PC,Mat*,Mat*,MatStructure*);
189 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
190 
191 extern PetscErrorCode PETSCKSP_DLLEXPORT PCView(PC,PetscViewer);
192 
193 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSetOptionsPrefix(PC,const char[]);
194 extern PetscErrorCode PETSCKSP_DLLEXPORT PCAppendOptionsPrefix(PC,const char[]);
195 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGetOptionsPrefix(PC,const char*[]);
196 
197 extern PetscErrorCode PETSCKSP_DLLEXPORT PCComputeExplicitOperator(PC,Mat*);
198 
199 /*
200       These are used to provide extra scaling of preconditioned
201    operator for time-stepping schemes like in SUNDIALS
202 */
203 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGetDiagonalScale(PC,PetscBool *);
204 extern PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleLeft(PC,Vec,Vec);
205 extern PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleRight(PC,Vec,Vec);
206 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSetDiagonalScale(PC,Vec);
207 
208 /* ------------- options specific to particular preconditioners --------- */
209 
210 extern PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowMax(PC);
211 extern PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowSum(PC);
212 extern PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseAbs(PC);
213 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric(PC,MatSORType);
214 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega(PC,PetscReal);
215 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations(PC,PetscInt,PetscInt);
216 
217 extern PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatSetOmega(PC,PetscReal);
218 extern PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatNoDiagonalScaling(PC);
219 
220 #define USE_PRECONDITIONER_MATRIX 0
221 #define USE_TRUE_MATRIX           1
222 extern PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetUseTrueLocal(PC);
223 extern PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
224 extern PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
225 
226 extern PetscErrorCode PETSCKSP_DLLEXPORT PCKSPSetUseTrue(PC);
227 
228 extern PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
229 extern PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
230 extern PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
231 extern PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
232 extern PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
233 extern PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
234 extern PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
235 extern PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetContext(PC,void**);
236 extern PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetContext(PC,void*);
237 extern PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetName(PC,const char[]);
238 extern PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetName(PC,char*[]);
239 
240 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot(PC,PetscReal);
241 
242 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftType(PC,MatFactorShiftType);
243 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftAmount(PC,PetscReal);
244 
245 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
246 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
247 
248 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC,PetscReal);
249 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetColumnPivot(PC,PetscReal);
250 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
251 
252 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC,const MatOrderingType);
253 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC,PetscBool );
254 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC,PetscBool );
255 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC);
256 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC);
257 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC,PetscBool );
258 
259 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC,PetscInt);
260 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
261 
262 extern PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
263 extern PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
264 extern PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetOverlap(PC,PetscInt);
265 extern PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetSortIndices(PC,PetscBool );
266 
267 /*E
268     PCASMType - Type of additive Schwarz method to use
269 
270 $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
271 $                 and computed values in ghost regions are added together. Classical
272 $                 standard additive Schwarz
273 $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
274 $                    region are discarded. Default
275 $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
276 $                       region are added back in
277 $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
278 $                not very good.
279 
280    Level: beginner
281 
282 .seealso: PCASMSetType()
283 E*/
284 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
285 extern const char *PCASMTypes[];
286 
287 extern PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetType(PC,PCASMType);
288 extern PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
289 extern PetscErrorCode PETSCKSP_DLLEXPORT PCASMDestroySubdomains(PetscInt,IS[],IS[]);
290 extern PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
291 extern PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
292 extern PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
293 
294 /*E
295     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per domain)
296 
297 $  PC_GASM_BASIC - symmetric version where residuals from the ghost points are used
298 $                 and computed values in ghost regions are added together. Classical
299 $                 standard additive Schwarz
300 $  PC_GASM_RESTRICT - residuals from ghost points are used but computed values in ghost
301 $                    region are discarded. Default
302 $  PC_GASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
303 $                       region are added back in
304 $  PC_GASM_NONE - ghost point residuals are not used, computed ghost values are discarded
305 $                not very good.
306 
307    Level: beginner
308 
309 .seealso: PCGASMSetType()
310 E*/
311 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
312 extern const char *PCGASMTypes[];
313 
314 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
315 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetTotalSubdomains(PC,PetscInt);
316 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetOverlap(PC,PetscInt);
317 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetSortIndices(PC,PetscBool );
318 
319 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGASMSetType(PC,PCGASMType);
320 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGASMCreateSubdomains(Mat,PetscInt,IS*[]);
321 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
322 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
323 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
324 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
325 
326 /*E
327     PCCompositeType - Determines how two or more preconditioner are composed
328 
329 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
330 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
331 $                                computed after the previous preconditioner application
332 $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
333 $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
334 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
335 $                         where first preconditioner is built from alpha I + S and second from
336 $                         alpha I + R
337 
338    Level: beginner
339 
340 .seealso: PCCompositeSetType()
341 E*/
342 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
343 extern const char *PCCompositeTypes[];
344 
345 extern PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC);
346 extern PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC,PCCompositeType);
347 extern PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC,PCType);
348 extern PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC,PetscInt,PC *);
349 extern PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC,PetscScalar);
350 
351 extern PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber(PC,PetscInt);
352 extern PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter(PC,VecScatter,VecScatter);
353 extern PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators(PC,Mat*,Mat*);
354 extern PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC(PC,PC*);
355 
356 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetEpsilon(PC,double);
357 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetNBSteps(PC,PetscInt);
358 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMax(PC,PetscInt);
359 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMaxNew(PC,PetscInt);
360 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetBlockSize(PC,PetscInt);
361 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetCacheSize(PC,PetscInt);
362 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetVerbose(PC,PetscInt);
363 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetSp(PC,PetscInt);
364 
365 extern PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC,const char[]);
366 extern PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC,const char*[]);
367 extern PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
368 extern PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
369 
370 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*);
371 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC,PCCompositeType);
372 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC,PetscInt);
373 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC,const char[],IS);
374 
375 /*E
376     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
377 
378     Level: intermediate
379 
380 .seealso: PCFieldSplitSchurPrecondition()
381 E*/
382 typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_DIAG,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
383 extern const char *const PCFieldSplitSchurPreTypes[];
384 
385 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
386 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
387 
388 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetRestriction(PC,Mat);
389 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetInterpolation(PC,Mat);
390 
391 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSetCoordinates(PC,PetscInt,PetscReal*);
392 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSASetVectors(PC,PetscInt,PetscReal *);
393 
394 extern PetscErrorCode PETSCKSP_DLLEXPORT PCPythonSetType(PC,const char[]);
395 
396 extern PetscErrorCode PETSCKSP_DLLEXPORT PCSetDM(PC,DM);
397 extern PetscErrorCode PETSCKSP_DLLEXPORT PCGetDM(PC,DM*);
398 
399 
400 PETSC_EXTERN_CXX_END
401 
402 #endif /* __PETSCPC_H */
403