xref: /petsc/include/petscpc.h (revision 3923b477fd0dced8a2d147b4fb4519fe3af97d3f)
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 /*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 PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
127 PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *);
128 PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool );
129 
130 PETSC_EXTERN PetscErrorCode PCRegisterDestroy(void);
131 PETSC_EXTERN PetscErrorCode PCRegisterAll(const char[]);
132 PETSC_EXTERN PetscBool PCRegisterAllCalled;
133 
134 PETSC_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 PETSC_EXTERN PetscErrorCode PCReset(PC);
184 PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
185 PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);
186 PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*);
187 
188 PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
189 PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
190 PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
191 
192 PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
193 PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
194 PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
195 
196 PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer);
197 
198 PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
199 PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
200 PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);
201 
202 PETSC_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 PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *);
209 PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
210 PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
211 PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec);
212 
213 /* ------------- options specific to particular preconditioners --------- */
214 
215 PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
216 PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC);
217 PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC);
218 PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
219 PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
220 PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
221 
222 PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
223 PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
224 
225 #define USE_PRECONDITIONER_MATRIX 0
226 #define USE_TRUE_MATRIX           1
227 PETSC_EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC);
228 PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
229 PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
230 
231 PETSC_EXTERN PetscErrorCode PCKSPSetUseTrue(PC);
232 
233 PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
234 PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
235 PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
236 PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
237 PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
238 PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
239 PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
240 PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**);
241 PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*);
242 PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
243 PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]);
244 
245 PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);
246 
247 PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType);
248 PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal);
249 
250 PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
251 PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
252 PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC);
253 
254 PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
255 PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal);
256 PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
257 
258 PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType);
259 PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool );
260 PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool );
261 PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC);
262 PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
263 PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool );
264 
265 PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
266 PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
267 
268 PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
269 PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
270 PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
271 PETSC_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.
278 $                        Classical standard additive Schwarz.
279 $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
280 $                        region are discarded.
281 $                        Default.
282 $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
283 $                        region are added back in.
284 $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
285 $                        discarded.
286 $                        Not very good.
287 
288    Level: beginner
289 
290 .seealso: PCASMSetType()
291 E*/
292 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
293 PETSC_EXTERN const char *const PCASMTypes[];
294 
295 PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
296 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
297 PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
298 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
299 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
300 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
301 
302 /*E
303     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
304 
305    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
306    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
307    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
308    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
309 
310 $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
311 $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
312 $                        from neighboring subdomains.
313 $                        Classical standard additive Schwarz.
314 $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
315 $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
316 $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
317 $                        assumption).
318 $                        Default.
319 $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
320 $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
321 $                        from neighboring subdomains.
322 $
323 $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
324 $                        Not very good.
325 
326    Level: beginner
327 
328 .seealso: PCGASMSetType()
329 E*/
330 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
331 PETSC_EXTERN const char *const PCGASMTypes[];
332 
333 PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
334 PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool);
335 PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
336 PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );
337 
338 PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
339 PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]);
340 PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
341 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
342 PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
343 PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);
344 
345 /*E
346     PCCompositeType - Determines how two or more preconditioner are composed
347 
348 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
349 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
350 $                                computed after the previous preconditioner application
351 $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
352 $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
353 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
354 $                         where first preconditioner is built from alpha I + S and second from
355 $                         alpha I + R
356 
357    Level: beginner
358 
359 .seealso: PCCompositeSetType()
360 E*/
361 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
362 PETSC_EXTERN const char *const PCCompositeTypes[];
363 
364 PETSC_EXTERN PetscErrorCode PCCompositeSetUseTrue(PC);
365 PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
366 PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
367 PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
368 PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
369 
370 PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
371 PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
372 PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
373 
374 PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
375 PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
376 PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
377 PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
378 PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
379 PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
380 PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
381 PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
382 
383 PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
384 PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
385 PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
386 PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
387 
388 PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
389 PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
390 PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
391 PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
392 PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
393 PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
394 
395 /*E
396     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
397 
398     Level: intermediate
399 
400 .seealso: PCFieldSplitSchurPrecondition()
401 E*/
402 typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_DIAG,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
403 PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
404 
405 /*E
406     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
407 
408     Level: intermediate
409 
410 .seealso: PCFieldSplitSetSchurFactType()
411 E*/
412 typedef enum {
413   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
414   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
415   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
416   PC_FIELDSPLIT_SCHUR_FACT_FULL
417 } PCFieldSplitSchurFactType;
418 PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
419 
420 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
421 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
422 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
423 
424 PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
425 PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
426 
427 PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
428 PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);
429 
430 PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);
431 
432 PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
433 PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);
434 
435 PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
436 PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);
437 
438 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
439 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
440 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
441 
442 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
443 PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
444 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
445 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
446 /*E
447     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
448 
449     Level: intermediate
450 
451 .seealso: PCPARMSSetGlobal()
452 E*/
453 typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
454 PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
455 /*E
456     PCPARMSLocalType - Determines the local preconditioner method in PARMS
457 
458     Level: intermediate
459 
460 .seealso: PCPARMSSetLocal()
461 E*/
462 typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
463 PETSC_EXTERN const char *const PCPARMSLocalTypes[];
464 
465 PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type);
466 PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type);
467 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits);
468 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart);
469 PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym);
470 PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2);
471 
472 PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
473 PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
474 PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
475 PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
476 PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
477 PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
478 PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
479 typedef const char* PCGAMGType;
480 PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType );
481 PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n);
482 PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n);
483 PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool);
484 
485 #if defined(PETSC_HAVE_PCBDDC)
486 /* Enum defining how to treat the coarse problem */
487 typedef enum {SEQUENTIAL_BDDC,REPLICATED_BDDC,PARALLEL_BDDC,MULTILEVEL_BDDC} CoarseProblemType;
488 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
489 PETSC_EXTERN PetscErrorCode PCBDDCSetMaxLevels(PC,PetscInt);
490 PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace);
491 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
492 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
493 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
494 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
495 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseProblemType(PC,CoarseProblemType);
496 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
497 PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,PetscInt[],PetscInt[],PetscCopyMode);
498 PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*);
499 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec);
500 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec);
501 #endif
502 
503 PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool);
504 PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
505 PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec);
506 
507 #endif /* __PETSCPC_H */
508