xref: /petsc/include/petscpc.h (revision 8e37d05fb63e2d28ce484af265873c94783d22f8)
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(void);
10 
11 /*
12     PCList contains the list of preconditioners currently registered
13    These are added with PCRegister()
14 */
15 PETSC_EXTERN PetscFunctionList PCList;
16 
17 /*S
18      PC - Abstract PETSc object that manages all preconditioners including direct solvers such as PCLU
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.
30 
31    Level: beginner
32 
33    Notes: Click on the links below to see details on a particular solver
34 
35           PCRegister() is used to register preconditioners that are then accessible via PCSetType()
36 
37 .seealso: PCSetType(), PC, PCCreate(), PCRegister(), PCSetFromOptions()
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 PCCP              "cp"
68 #define PCBFBT            "bfbt"
69 #define PCLSC             "lsc"
70 #define PCPYTHON          "python"
71 #define PCPFMG            "pfmg"
72 #define PCSYSPFMG         "syspfmg"
73 #define PCREDISTRIBUTE    "redistribute"
74 #define PCSVD             "svd"
75 #define PCGAMG            "gamg"
76 #define PCSACUSP          "sacusp"        /* these four run on NVIDIA GPUs using CUSP */
77 #define PCSACUSPPOLY      "sacusppoly"
78 #define PCBICGSTABCUSP    "bicgstabcusp"
79 #define PCAINVCUSP        "ainvcusp"
80 #define PCBDDC            "bddc"
81 #define PCKACZMARZ        "kaczmarz"
82 
83 /* Logging support */
84 PETSC_EXTERN PetscClassId PC_CLASSID;
85 
86 /*E
87     PCSide - If the preconditioner is to be applied to the left, right
88      or symmetrically around the operator.
89 
90    Level: beginner
91 
92 .seealso:
93 E*/
94 typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
95 #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
96 PETSC_EXTERN const char *const *const PCSides;
97 
98 PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
99 PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType);
100 PETSC_EXTERN PetscErrorCode PCSetUp(PC);
101 PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
102 PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
103 PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
104 PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
105 PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
106 PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
107 PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *);
108 PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
109 PETSC_EXTERN PetscErrorCode PCSetReusePreconditioner(PC,PetscBool);
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 PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool);
132 PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*);
133 
134 PETSC_EXTERN PetscErrorCode PCRegisterAll(void);
135 PETSC_EXTERN PetscBool PCRegisterAllCalled;
136 
137 PETSC_EXTERN PetscErrorCode PCRegister(const char[],PetscErrorCode(*)(PC));
138 
139 PETSC_EXTERN PetscErrorCode PCReset(PC);
140 PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
141 PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);
142 PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*);
143 
144 PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
145 PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
146 PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
147 
148 PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat);
149 PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*);
150 PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
151 
152 PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer);
153 PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer);
154 PETSC_STATIC_INLINE PetscErrorCode PCViewFromOptions(PC A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
155 
156 PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
157 PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
158 PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);
159 
160 PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
161 
162 /*
163       These are used to provide extra scaling of preconditioned
164    operator for time-stepping schemes like in SUNDIALS
165 */
166 PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *);
167 PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
168 PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
169 PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec);
170 
171 /* ------------- options specific to particular preconditioners --------- */
172 
173 /*E
174     PCJacobiType - What elements are used to form the Jacobi preconditioner
175 
176    Level: intermediate
177 
178 .seealso:
179 E*/
180 typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType;
181 PETSC_EXTERN const char *const PCJacobiTypes[];
182 
183 PETSC_EXTERN PetscErrorCode PCJacobiSetType(PC,PCJacobiType);
184 PETSC_EXTERN PetscErrorCode PCJacobiGetType(PC,PCJacobiType*);
185 PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC,PetscBool);
186 PETSC_EXTERN PetscErrorCode PCJacobiGetUseAbs(PC,PetscBool*);
187 PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
188 PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
189 PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
190 
191 PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
192 PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
193 
194 PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
195 PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
196 
197 PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
198 PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
199 PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
200 PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
201 PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
202 PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
203 PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
204 PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**);
205 PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*);
206 PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
207 PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]);
208 
209 PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);
210 
211 PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType);
212 PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal);
213 
214 PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
215 PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
216 PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC);
217 
218 PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
219 PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal);
220 PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
221 
222 PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType);
223 PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool );
224 PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool );
225 PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC,PetscBool);
226 PETSC_EXTERN PetscErrorCode PCFactorGetUseInPlace(PC,PetscBool*);
227 PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
228 PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool );
229 
230 PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC,PetscInt*);
231 PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
232 PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
233 
234 PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
235 PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
236 PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
237 PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool);
238 PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*);
239 PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool );
240 
241 /*E
242     PCASMType - Type of additive Schwarz method to use
243 
244 $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
245 $                        and computed values in ghost regions are added together.
246 $                        Classical standard additive Schwarz.
247 $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
248 $                        region are discarded.
249 $                        Default.
250 $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
251 $                        region are added back in.
252 $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
253 $                        discarded.
254 $                        Not very good.
255 
256    Level: beginner
257 
258 .seealso: PCASMSetType()
259 E*/
260 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
261 PETSC_EXTERN const char *const PCASMTypes[];
262 
263 PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
264 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
265 PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
266 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
267 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
268 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
269 
270 /*E
271     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
272 
273    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
274    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
275    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
276    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
277 
278 $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
279 $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
280 $                        from neighboring subdomains.
281 $                        Classical standard additive Schwarz.
282 $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
283 $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
284 $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
285 $                        assumption).
286 $                        Default.
287 $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
288 $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
289 $                        from neighboring subdomains.
290 $
291 $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
292 $                        Not very good.
293 
294    Level: beginner
295 
296 .seealso: PCGASMSetType()
297 E*/
298 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
299 PETSC_EXTERN const char *const PCGASMTypes[];
300 
301 PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
302 PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool);
303 PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
304 PETSC_EXTERN PetscErrorCode PCGASMSetDMSubdomains(PC,PetscBool);
305 PETSC_EXTERN PetscErrorCode PCGASMGetDMSubdomains(PC,PetscBool*);
306 PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );
307 
308 PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
309 PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]);
310 PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
311 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
312 PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
313 PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);
314 
315 /*E
316     PCCompositeType - Determines how two or more preconditioner are composed
317 
318 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
319 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
320 $                                computed after the previous preconditioner application
321 $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
322 $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
323 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
324 $                         where first preconditioner is built from alpha I + S and second from
325 $                         alpha I + R
326 
327    Level: beginner
328 
329 .seealso: PCCompositeSetType()
330 E*/
331 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
332 PETSC_EXTERN const char *const PCCompositeTypes[];
333 
334 PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
335 PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
336 PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
337 PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
338 
339 PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
340 PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
341 PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
342 
343 PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
344 PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
345 PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
346 PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
347 PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
348 PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
349 PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
350 PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
351 
352 PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
353 PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
354 PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
355 PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
356 
357 PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
358 PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
359 PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
360 PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
361 PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
362 PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
363 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool);
364 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*);
365 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC,PetscBool);
366 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC,PetscBool*);
367 PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC,PetscBool);
368 PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC,PetscBool*);
369 
370 
371 /*E
372     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
373 
374     Level: intermediate
375 
376 .seealso: PCFieldSplitSetSchurPre()
377 E*/
378 typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_SELFP,PC_FIELDSPLIT_SCHUR_PRE_A11,PC_FIELDSPLIT_SCHUR_PRE_USER,PC_FIELDSPLIT_SCHUR_PRE_FULL} PCFieldSplitSchurPreType;
379 PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
380 
381 /*E
382     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
383 
384     Level: intermediate
385 
386 .seealso: PCFieldSplitSetSchurFactType()
387 E*/
388 typedef enum {
389   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
390   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
391   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
392   PC_FIELDSPLIT_SCHUR_FACT_FULL
393 } PCFieldSplitSchurFactType;
394 PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
395 
396 PETSC_EXTERN PETSC_DEPRECATED("Use PCFieldSplitSetSchurPre") PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
397 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC,PCFieldSplitSchurPreType,Mat);
398 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC,PCFieldSplitSchurPreType*,Mat*);
399 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
400 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
401 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC,Mat *S);
402 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC,Mat *S);
403 
404 PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
405 PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
406 
407 PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
408 
409 PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);
410 
411 PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
412 PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);
413 
414 PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
415 PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);
416 
417 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
418 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
419 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
420 
421 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
422 PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
423 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
424 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
425 /*E
426     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
427 
428     Level: intermediate
429 
430 .seealso: PCPARMSSetGlobal()
431 E*/
432 typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
433 PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
434 /*E
435     PCPARMSLocalType - Determines the local preconditioner method in PARMS
436 
437     Level: intermediate
438 
439 .seealso: PCPARMSSetLocal()
440 E*/
441 typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
442 PETSC_EXTERN const char *const PCPARMSLocalTypes[];
443 
444 PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC,PCPARMSGlobalType);
445 PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC,PCPARMSLocalType);
446 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC,PetscReal,PetscInt);
447 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC,PetscInt);
448 PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC,PetscBool);
449 PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC,PetscInt,PetscInt,PetscInt);
450 
451 /*E
452     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method
453 
454     Level: intermediate
455 
456 .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseProl()
457 E*/
458 typedef const char *PCGAMGType;
459 #define PCGAMGAGG         "agg"
460 #define PCGAMGGEO         "geo"
461 #define PCGAMGCLASSICAL   "classical"
462 PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
463 PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
464 PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
465 PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
466 PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
467 PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
468 PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
469 PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType);
470 PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC,PetscInt);
471 PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC,PetscBool);
472 PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool);
473 PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool);
474 PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void);
475 PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void);
476 
477 typedef const char *PCGAMGClassicalType;
478 #define PCGAMGCLASSICALDIRECT   "direct"
479 #define PCGAMGCLASSICALSTANDARD "standard"
480 PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC,PCGAMGClassicalType);
481 
482 #if defined(PETSC_HAVE_PCBDDC)
483 PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisMat(PC,Mat);
484 PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC,IS);
485 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
486 PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC,PetscInt);
487 PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace);
488 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
489 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC,IS);
490 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
491 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC,IS*);
492 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
493 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC,IS);
494 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
495 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC,IS*);
496 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
497 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC,PetscInt,IS[]);
498 PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode);
499 PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*);
500 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec);
501 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec);
502 #endif
503 
504 PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool);
505 PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
506 PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec);
507 
508 /*E
509     PCMGType - Determines the type of multigrid method that is run.
510 
511    Level: beginner
512 
513    Values:
514 +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles()
515 .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
516                 smoothed before updating the residual. This only uses the
517                 down smoother, in the preconditioner the upper smoother is ignored
518 .  PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
519             that is starts on the coarsest grid, performs a cycle, interpolates
520             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
521             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
522             calls the V-cycle only on the coarser level and has a post-smoother instead.
523 -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
524                from a finer
525 
526 .seealso: PCMGSetType()
527 
528 E*/
529 typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
530 PETSC_EXTERN const char *const PCMGTypes[];
531 #define PC_MG_CASCADE PC_MG_KASKADE;
532 
533 /*E
534     PCMGCycleType - Use V-cycle or W-cycle
535 
536    Level: beginner
537 
538    Values:
539 +  PC_MG_V_CYCLE
540 -  PC_MG_W_CYCLE
541 
542 .seealso: PCMGSetCycleType()
543 
544 E*/
545 typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
546 PETSC_EXTERN const char *const PCMGCycleTypes[];
547 
548 PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType);
549 PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*);
550 PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*);
551 
552 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt);
553 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt);
554 PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType);
555 PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType);
556 PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt);
557 PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt);
558 PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool);
559 PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*);
560 
561 
562 PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec);
563 PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec);
564 PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec);
565 
566 PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat);
567 PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*);
568 PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat);
569 PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*);
570 PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec);
571 PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*);
572 PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
573 PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat,Vec,Vec,Vec);
574 
575 /*E
576     PCExoticType - Face based or wirebasket based coarse grid space
577 
578    Level: beginner
579 
580 .seealso: PCExoticSetType(), PCEXOTIC
581 E*/
582 typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
583 PETSC_EXTERN const char *const PCExoticTypes[];
584 PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);
585 
586 #endif /* __PETSCPC_H */
587