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