xref: /petsc/include/petscpc.h (revision 242ebda14a69a344bed43bec2e001e417d7cf4e9)
1 /*
2       Preconditioner module.
3 */
4 #if !defined(__PETSCPC_H)
5 #define __PETSCPC_H
6 #include <petscmat.h>
7 #include <petscdmtypes.h>
8 
9 PETSC_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 PETSC_EXTERN PetscFunctionList 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 /*J
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 J*/
39 typedef const char* PCType;
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 PCGALERKIN        "galerkin"
66 #define PCEXOTIC          "exotic"
67 #define PCHMPI            "hmpi"
68 #define PCSUPPORTGRAPH    "supportgraph"
69 #define PCASA             "asa"
70 #define PCCP              "cp"
71 #define PCBFBT            "bfbt"
72 #define PCLSC             "lsc"
73 #define PCPYTHON          "python"
74 #define PCPFMG            "pfmg"
75 #define PCSYSPFMG         "syspfmg"
76 #define PCREDISTRIBUTE    "redistribute"
77 #define PCSVD             "svd"
78 #define PCGAMG            "gamg"
79 #define PCSACUSP          "sacusp"        /* these four run on NVIDIA GPUs using CUSP */
80 #define PCSACUSPPOLY      "sacusppoly"
81 #define PCBICGSTABCUSP    "bicgstabcusp"
82 #define PCAINVCUSP        "ainvcusp"
83 #define PCBDDC            "bddc"
84 
85 /* Logging support */
86 PETSC_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_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
97 #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
98 PETSC_EXTERN const char *const *const PCSides;
99 
100 PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
101 PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType);
102 PETSC_EXTERN PetscErrorCode PCSetUp(PC);
103 PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
104 PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
105 PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
106 PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
107 PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
108 PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
109 PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *);
110 PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
111 
112 #define PC_FILE_CLASSID 1211222
113 
114 /*E
115     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
116 
117    Level: advanced
118 
119    Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
120 
121 .seealso: PCApplyRichardson()
122 E*/
123 typedef enum {
124               PCRICHARDSON_CONVERGED_RTOL               =  2,
125               PCRICHARDSON_CONVERGED_ATOL               =  3,
126               PCRICHARDSON_CONVERGED_ITS                =  4,
127               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
128 
129 PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
130 PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *);
131 PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool );
132 PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool);
133 PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*);
134 
135 PETSC_EXTERN PetscErrorCode PCRegisterDestroy(void);
136 PETSC_EXTERN PetscErrorCode PCRegisterAll(const char[]);
137 PETSC_EXTERN PetscBool PCRegisterAllCalled;
138 
139 PETSC_EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
140 
141 /*MC
142    PCRegisterDynamic - Adds a method to the preconditioner package.
143 
144    Synopsis:
145     #include "petscpc.h"
146    PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC))
147 
148    Not collective
149 
150    Input Parameters:
151 +  name_solver - name of a new user-defined solver
152 .  path - path (either absolute or relative) the library containing this solver
153 .  name_create - name of routine to create method context
154 -  routine_create - routine to create method context
155 
156    Notes:
157    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
158 
159    If dynamic libraries are used, then the fourth input argument (routine_create)
160    is ignored.
161 
162    Sample usage:
163 .vb
164    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
165               "MySolverCreate",MySolverCreate);
166 .ve
167 
168    Then, your solver can be chosen with the procedural interface via
169 $     PCSetType(pc,"my_solver")
170    or at runtime via the option
171 $     -pc_type my_solver
172 
173    Level: advanced
174 
175    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
176            occuring in pathname will be replaced with appropriate values.
177          If your function is not being put into a shared library then use PCRegister() instead
178 
179 .keywords: PC, register
180 
181 .seealso: PCRegisterAll(), PCRegisterDestroy()
182 M*/
183 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
184 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
185 #else
186 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
187 #endif
188 
189 PETSC_EXTERN PetscErrorCode PCReset(PC);
190 PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
191 PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);
192 PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*);
193 
194 PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
195 PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
196 PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
197 
198 PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
199 PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
200 PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
201 
202 PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer);
203 PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer);
204 
205 PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
206 PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
207 PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);
208 
209 PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
210 
211 /*
212       These are used to provide extra scaling of preconditioned
213    operator for time-stepping schemes like in SUNDIALS
214 */
215 PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *);
216 PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
217 PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
218 PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec);
219 
220 /* ------------- options specific to particular preconditioners --------- */
221 
222 PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
223 PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC);
224 PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC);
225 PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
226 PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
227 PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
228 
229 PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
230 PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
231 
232 PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
233 PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
234 
235 PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
236 PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
237 PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
238 PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
239 PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
240 PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
241 PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
242 PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**);
243 PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*);
244 PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
245 PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]);
246 
247 PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);
248 
249 PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType);
250 PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal);
251 
252 PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
253 PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
254 PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC);
255 
256 PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
257 PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal);
258 PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
259 
260 PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType);
261 PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool );
262 PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool );
263 PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC);
264 PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
265 PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool );
266 
267 PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
268 PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
269 
270 PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
271 PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
272 PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
273 PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool);
274 PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*);
275 PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool );
276 
277 /*E
278     PCASMType - Type of additive Schwarz method to use
279 
280 $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
281 $                        and computed values in ghost regions are added together.
282 $                        Classical standard additive Schwarz.
283 $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
284 $                        region are discarded.
285 $                        Default.
286 $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
287 $                        region are added back in.
288 $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
289 $                        discarded.
290 $                        Not very good.
291 
292    Level: beginner
293 
294 .seealso: PCASMSetType()
295 E*/
296 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
297 PETSC_EXTERN const char *const PCASMTypes[];
298 
299 PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
300 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
301 PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
302 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
303 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
304 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
305 
306 /*E
307     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
308 
309    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
310    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
311    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
312    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
313 
314 $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
315 $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
316 $                        from neighboring subdomains.
317 $                        Classical standard additive Schwarz.
318 $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
319 $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
320 $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
321 $                        assumption).
322 $                        Default.
323 $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
324 $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
325 $                        from neighboring subdomains.
326 $
327 $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
328 $                        Not very good.
329 
330    Level: beginner
331 
332 .seealso: PCGASMSetType()
333 E*/
334 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
335 PETSC_EXTERN const char *const PCGASMTypes[];
336 
337 PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
338 PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool);
339 PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
340 PETSC_EXTERN PetscErrorCode PCGASMSetDMSubdomains(PC,PetscBool);
341 PETSC_EXTERN PetscErrorCode PCGASMGetDMSubdomains(PC,PetscBool*);
342 PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );
343 
344 PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
345 PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]);
346 PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
347 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
348 PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
349 PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);
350 
351 /*E
352     PCCompositeType - Determines how two or more preconditioner are composed
353 
354 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
355 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
356 $                                computed after the previous preconditioner application
357 $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
358 $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
359 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
360 $                         where first preconditioner is built from alpha I + S and second from
361 $                         alpha I + R
362 
363    Level: beginner
364 
365 .seealso: PCCompositeSetType()
366 E*/
367 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
368 PETSC_EXTERN const char *const PCCompositeTypes[];
369 
370 PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
371 PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
372 PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
373 PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
374 
375 PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
376 PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
377 PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
378 
379 PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
380 PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
381 PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
382 PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
383 PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
384 PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
385 PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
386 PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
387 
388 PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
389 PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
390 PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
391 PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
392 
393 PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
394 PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
395 PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
396 PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
397 PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
398 PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
399 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool);
400 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*);
401 
402 /*E
403     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
404 
405     Level: intermediate
406 
407 .seealso: PCFieldSplitSchurPrecondition()
408 E*/
409 typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_A11,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
410 PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
411 
412 /*E
413     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
414 
415     Level: intermediate
416 
417 .seealso: PCFieldSplitSetSchurFactType()
418 E*/
419 typedef enum {
420   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
421   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
422   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
423   PC_FIELDSPLIT_SCHUR_FACT_FULL
424 } PCFieldSplitSchurFactType;
425 PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
426 
427 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
428 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
429 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
430 
431 PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
432 PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
433 
434 PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
435 PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);
436 
437 PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);
438 
439 PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
440 PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);
441 
442 PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
443 PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);
444 
445 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
446 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
447 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
448 
449 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
450 PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
451 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
452 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
453 /*E
454     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
455 
456     Level: intermediate
457 
458 .seealso: PCPARMSSetGlobal()
459 E*/
460 typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
461 PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
462 /*E
463     PCPARMSLocalType - Determines the local preconditioner method in PARMS
464 
465     Level: intermediate
466 
467 .seealso: PCPARMSSetLocal()
468 E*/
469 typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
470 PETSC_EXTERN const char *const PCPARMSLocalTypes[];
471 
472 PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type);
473 PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type);
474 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits);
475 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart);
476 PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym);
477 PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2);
478 
479 PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
480 PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
481 PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
482 PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
483 PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
484 PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
485 PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
486 typedef const char* PCGAMGType;
487 PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType );
488 PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n);
489 PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n);
490 PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool);
491 PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool);
492 
493 #if defined(PETSC_HAVE_PCBDDC)
494 /* Enum defining how to treat the coarse problem */
495 typedef enum {SEQUENTIAL_BDDC,REPLICATED_BDDC,PARALLEL_BDDC,MULTILEVEL_BDDC} CoarseProblemType;
496 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
497 PETSC_EXTERN PetscErrorCode PCBDDCSetMaxLevels(PC,PetscInt);
498 PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace);
499 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
500 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
501 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
502 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
503 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseProblemType(PC,CoarseProblemType);
504 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
505 PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode);
506 PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*);
507 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec);
508 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec);
509 #endif
510 
511 PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool);
512 PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
513 PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec);
514 
515 /*E
516     PCMGType - Determines the type of multigrid method that is run.
517 
518    Level: beginner
519 
520    Values:
521 +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles()
522 .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
523                 smoothed before updating the residual. This only uses the
524                 down smoother, in the preconditioner the upper smoother is ignored
525 .  PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
526             that is starts on the coarsest grid, performs a cycle, interpolates
527             to the next, performs a cycle etc. This is much like the F-cycle presented in "Multigrid" by Trottenberg, Oosterlee, Schuller page 49, but that
528             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
529             calls the V-cycle only on the coarser level and has a post-smoother instead.
530 -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
531                from a finer
532 
533 .seealso: PCMGSetType()
534 
535 E*/
536 typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
537 PETSC_EXTERN const char *const PCMGTypes[];
538 #define PC_MG_CASCADE PC_MG_KASKADE;
539 
540 /*E
541     PCMGCycleType - Use V-cycle or W-cycle
542 
543    Level: beginner
544 
545    Values:
546 +  PC_MG_V_CYCLE
547 -  PC_MG_W_CYCLE
548 
549 .seealso: PCMGSetCycleType()
550 
551 E*/
552 typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
553 PETSC_EXTERN const char *const PCMGCycleTypes[];
554 
555 PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType);
556 PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*);
557 PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*);
558 
559 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt);
560 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt);
561 PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType);
562 PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType);
563 PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt);
564 PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt);
565 PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool);
566 PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*);
567 
568 
569 PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec);
570 PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec);
571 PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec);
572 
573 PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat);
574 PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*);
575 PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat);
576 PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*);
577 PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec);
578 PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*);
579 PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
580 
581 /*E
582     PCExoticType - Face based or wirebasket based coarse grid space
583 
584    Level: beginner
585 
586 .seealso: PCExoticSetType(), PCEXOTIC
587 E*/
588 typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
589 PETSC_EXTERN const char *const PCExoticTypes[];
590 PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);
591 
592 #endif /* __PETSCPC_H */
593