xref: /petsc/include/petscpc.h (revision f2c2a1b97ba623b07d38ea716c110451b7b5ea3c)
1 /*
2       Preconditioner module.
3 */
4 #if !defined(__PETSCPC_H)
5 #define __PETSCPC_H
6 #include <petscdm.h>
7 
8 PETSC_EXTERN PetscErrorCode PCInitializePackage(const char[]);
9 
10 /*
11     PCList contains the list of preconditioners currently registered
12    These are added with the PCRegisterDynamic() macro
13 */
14 PETSC_EXTERN PetscFList PCList;
15 
16 /*S
17      PC - Abstract PETSc object that manages all preconditioners
18 
19    Level: beginner
20 
21   Concepts: preconditioners
22 
23 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
24 S*/
25 typedef struct _p_PC* PC;
26 
27 /*J
28     PCType - String with the name of a PETSc preconditioner method or the creation function
29        with an optional dynamic library name, for example
30        http://www.mcs.anl.gov/petsc/lib.a:mypccreate()
31 
32    Level: beginner
33 
34    Notes: Click on the links below to see details on a particular solver
35 
36 .seealso: PCSetType(), PC, PCCreate()
37 J*/
38 typedef const char* PCType;
39 #define PCNONE            "none"
40 #define PCJACOBI          "jacobi"
41 #define PCSOR             "sor"
42 #define PCLU              "lu"
43 #define PCSHELL           "shell"
44 #define PCBJACOBI         "bjacobi"
45 #define PCMG              "mg"
46 #define PCEISENSTAT       "eisenstat"
47 #define PCILU             "ilu"
48 #define PCICC             "icc"
49 #define PCASM             "asm"
50 #define PCGASM            "gasm"
51 #define PCKSP             "ksp"
52 #define PCCOMPOSITE       "composite"
53 #define PCREDUNDANT       "redundant"
54 #define PCSPAI            "spai"
55 #define PCNN              "nn"
56 #define PCCHOLESKY        "cholesky"
57 #define PCPBJACOBI        "pbjacobi"
58 #define PCMAT             "mat"
59 #define PCHYPRE           "hypre"
60 #define PCPARMS           "parms"
61 #define PCFIELDSPLIT      "fieldsplit"
62 #define PCTFS             "tfs"
63 #define PCML              "ml"
64 #define PCGALERKIN        "galerkin"
65 #define PCEXOTIC          "exotic"
66 #define PCHMPI            "hmpi"
67 #define PCSUPPORTGRAPH    "supportgraph"
68 #define PCASA             "asa"
69 #define PCCP              "cp"
70 #define PCBFBT            "bfbt"
71 #define PCLSC             "lsc"
72 #define PCPYTHON          "python"
73 #define PCPFMG            "pfmg"
74 #define PCSYSPFMG         "syspfmg"
75 #define PCREDISTRIBUTE    "redistribute"
76 #define PCSVD             "svd"
77 #define PCGAMG            "gamg"
78 #define PCSACUSP          "sacusp"        /* these four run on NVIDIA GPUs using CUSP */
79 #define PCSACUSPPOLY      "sacusppoly"
80 #define PCBICGSTABCUSP    "bicgstabcusp"
81 #define PCAINVCUSP        "ainvcusp"
82 #define PCBDDC            "bddc"
83 
84 /* Logging support */
85 PETSC_EXTERN PetscClassId PC_CLASSID;
86 
87 /*E
88     PCSide - If the preconditioner is to be applied to the left, right
89      or symmetrically around the operator.
90 
91    Level: beginner
92 
93 .seealso:
94 E*/
95 typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
96 #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
97 PETSC_EXTERN const char *const *const PCSides;
98 
99 PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
100 PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType);
101 PETSC_EXTERN PetscErrorCode PCSetUp(PC);
102 PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
103 PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
104 PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
105 PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
106 PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
107 PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
108 PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *);
109 PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
110 
111 #define PC_FILE_CLASSID 1211222
112 
113 /*E
114     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
115 
116    Level: advanced
117 
118    Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
119 
120 .seealso: PCApplyRichardson()
121 E*/
122 typedef enum {
123               PCRICHARDSON_CONVERGED_RTOL               =  2,
124               PCRICHARDSON_CONVERGED_ATOL               =  3,
125               PCRICHARDSON_CONVERGED_ITS                =  4,
126               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
127 
128 PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
129 PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *);
130 PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool );
131 
132 PETSC_EXTERN PetscErrorCode PCRegisterDestroy(void);
133 PETSC_EXTERN PetscErrorCode PCRegisterAll(const char[]);
134 PETSC_EXTERN PetscBool PCRegisterAllCalled;
135 
136 PETSC_EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
137 
138 /*MC
139    PCRegisterDynamic - Adds a method to the preconditioner package.
140 
141    Synopsis:
142    PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC))
143 
144    Not collective
145 
146    Input Parameters:
147 +  name_solver - name of a new user-defined solver
148 .  path - path (either absolute or relative) the library containing this solver
149 .  name_create - name of routine to create method context
150 -  routine_create - routine to create method context
151 
152    Notes:
153    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
154 
155    If dynamic libraries are used, then the fourth input argument (routine_create)
156    is ignored.
157 
158    Sample usage:
159 .vb
160    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
161               "MySolverCreate",MySolverCreate);
162 .ve
163 
164    Then, your solver can be chosen with the procedural interface via
165 $     PCSetType(pc,"my_solver")
166    or at runtime via the option
167 $     -pc_type my_solver
168 
169    Level: advanced
170 
171    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
172            occuring in pathname will be replaced with appropriate values.
173          If your function is not being put into a shared library then use PCRegister() instead
174 
175 .keywords: PC, register
176 
177 .seealso: PCRegisterAll(), PCRegisterDestroy()
178 M*/
179 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
180 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
181 #else
182 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
183 #endif
184 
185 PETSC_EXTERN PetscErrorCode PCReset(PC);
186 PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
187 PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);
188 PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*);
189 
190 PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
191 PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
192 PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
193 
194 PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
195 PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
196 PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
197 
198 PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer);
199 PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer);
200 
201 PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
202 PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
203 PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);
204 
205 PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
206 
207 /*
208       These are used to provide extra scaling of preconditioned
209    operator for time-stepping schemes like in SUNDIALS
210 */
211 PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *);
212 PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
213 PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
214 PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec);
215 
216 /* ------------- options specific to particular preconditioners --------- */
217 
218 PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
219 PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC);
220 PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC);
221 PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
222 PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
223 PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
224 
225 PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
226 PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
227 
228 #define USE_PRECONDITIONER_MATRIX 0
229 #define USE_TRUE_MATRIX           1
230 PETSC_EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC);
231 PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
232 PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
233 
234 PETSC_EXTERN PetscErrorCode PCKSPSetUseTrue(PC);
235 
236 PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
237 PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
238 PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
239 PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
240 PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
241 PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
242 PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
243 PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**);
244 PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*);
245 PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
246 PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]);
247 
248 PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);
249 
250 PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType);
251 PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal);
252 
253 PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
254 PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
255 PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC);
256 
257 PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
258 PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal);
259 PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
260 
261 PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType);
262 PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool );
263 PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool );
264 PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC);
265 PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
266 PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool );
267 
268 PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
269 PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
270 
271 PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
272 PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
273 PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
274 PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool );
275 
276 /*E
277     PCASMType - Type of additive Schwarz method to use
278 
279 $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
280 $                        and computed values in ghost regions are added together.
281 $                        Classical standard additive Schwarz.
282 $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
283 $                        region are discarded.
284 $                        Default.
285 $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
286 $                        region are added back in.
287 $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
288 $                        discarded.
289 $                        Not very good.
290 
291    Level: beginner
292 
293 .seealso: PCASMSetType()
294 E*/
295 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
296 PETSC_EXTERN const char *const PCASMTypes[];
297 
298 PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
299 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
300 PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
301 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
302 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
303 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
304 
305 /*E
306     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
307 
308    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
309    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
310    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
311    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
312 
313 $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
314 $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
315 $                        from neighboring subdomains.
316 $                        Classical standard additive Schwarz.
317 $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
318 $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
319 $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
320 $                        assumption).
321 $                        Default.
322 $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
323 $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
324 $                        from neighboring subdomains.
325 $
326 $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
327 $                        Not very good.
328 
329    Level: beginner
330 
331 .seealso: PCGASMSetType()
332 E*/
333 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
334 PETSC_EXTERN const char *const PCGASMTypes[];
335 
336 PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
337 PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool);
338 PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
339 PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );
340 
341 PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
342 PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]);
343 PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
344 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
345 PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
346 PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);
347 
348 /*E
349     PCCompositeType - Determines how two or more preconditioner are composed
350 
351 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
352 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
353 $                                computed after the previous preconditioner application
354 $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
355 $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
356 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
357 $                         where first preconditioner is built from alpha I + S and second from
358 $                         alpha I + R
359 
360    Level: beginner
361 
362 .seealso: PCCompositeSetType()
363 E*/
364 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
365 PETSC_EXTERN const char *const PCCompositeTypes[];
366 
367 PETSC_EXTERN PetscErrorCode PCCompositeSetUseTrue(PC);
368 PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
369 PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
370 PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
371 PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
372 
373 PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
374 PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
375 PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
376 
377 PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
378 PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
379 PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
380 PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
381 PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
382 PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
383 PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
384 PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
385 
386 PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
387 PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
388 PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
389 PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
390 
391 PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
392 PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
393 PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
394 PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
395 PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
396 PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
397 
398 /*E
399     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
400 
401     Level: intermediate
402 
403 .seealso: PCFieldSplitSchurPrecondition()
404 E*/
405 typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_DIAG,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
406 PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
407 
408 /*E
409     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
410 
411     Level: intermediate
412 
413 .seealso: PCFieldSplitSetSchurFactType()
414 E*/
415 typedef enum {
416   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
417   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
418   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
419   PC_FIELDSPLIT_SCHUR_FACT_FULL
420 } PCFieldSplitSchurFactType;
421 PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
422 
423 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
424 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
425 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
426 
427 PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
428 PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
429 
430 PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
431 PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);
432 
433 PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);
434 
435 PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
436 PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);
437 
438 PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
439 PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);
440 
441 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
442 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
443 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
444 
445 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
446 PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
447 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
448 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
449 /*E
450     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
451 
452     Level: intermediate
453 
454 .seealso: PCPARMSSetGlobal()
455 E*/
456 typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
457 PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
458 /*E
459     PCPARMSLocalType - Determines the local preconditioner method in PARMS
460 
461     Level: intermediate
462 
463 .seealso: PCPARMSSetLocal()
464 E*/
465 typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
466 PETSC_EXTERN const char *const PCPARMSLocalTypes[];
467 
468 PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type);
469 PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type);
470 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits);
471 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart);
472 PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym);
473 PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2);
474 
475 PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
476 PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
477 PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
478 PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
479 PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
480 PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
481 PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
482 typedef const char* PCGAMGType;
483 PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType );
484 PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n);
485 PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n);
486 PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool);
487 
488 #if defined(PETSC_HAVE_PCBDDC)
489 /* Enum defining how to treat the coarse problem */
490 typedef enum {SEQUENTIAL_BDDC,REPLICATED_BDDC,PARALLEL_BDDC,MULTILEVEL_BDDC} CoarseProblemType;
491 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
492 PETSC_EXTERN PetscErrorCode PCBDDCSetMaxLevels(PC,PetscInt);
493 PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace);
494 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
495 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
496 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
497 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
498 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseProblemType(PC,CoarseProblemType);
499 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
500 PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode);
501 PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*);
502 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec);
503 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec);
504 #endif
505 
506 PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool);
507 PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
508 PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec);
509 
510 #endif /* __PETSCPC_H */
511