xref: /petsc/include/petscpc.h (revision baa89ecb115416e2a197e86cd1cb4e9abb543814)
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);
226 PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
227 PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool );
228 
229 PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC,PetscInt*);
230 PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
231 PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
232 
233 PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
234 PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
235 PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
236 PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool);
237 PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*);
238 PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool );
239 
240 /*E
241     PCASMType - Type of additive Schwarz method to use
242 
243 $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
244 $                        and computed values in ghost regions are added together.
245 $                        Classical standard additive Schwarz.
246 $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
247 $                        region are discarded.
248 $                        Default.
249 $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
250 $                        region are added back in.
251 $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
252 $                        discarded.
253 $                        Not very good.
254 
255    Level: beginner
256 
257 .seealso: PCASMSetType()
258 E*/
259 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
260 PETSC_EXTERN const char *const PCASMTypes[];
261 
262 PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
263 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
264 PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
265 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
266 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
267 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
268 
269 /*E
270     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
271 
272    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
273    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
274    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
275    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
276 
277 $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
278 $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
279 $                        from neighboring subdomains.
280 $                        Classical standard additive Schwarz.
281 $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
282 $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
283 $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
284 $                        assumption).
285 $                        Default.
286 $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
287 $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
288 $                        from neighboring subdomains.
289 $
290 $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
291 $                        Not very good.
292 
293    Level: beginner
294 
295 .seealso: PCGASMSetType()
296 E*/
297 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
298 PETSC_EXTERN const char *const PCGASMTypes[];
299 
300 PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
301 PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool);
302 PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
303 PETSC_EXTERN PetscErrorCode PCGASMSetDMSubdomains(PC,PetscBool);
304 PETSC_EXTERN PetscErrorCode PCGASMGetDMSubdomains(PC,PetscBool*);
305 PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );
306 
307 PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
308 PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]);
309 PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
310 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
311 PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
312 PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);
313 
314 /*E
315     PCCompositeType - Determines how two or more preconditioner are composed
316 
317 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
318 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
319 $                                computed after the previous preconditioner application
320 $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
321 $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
322 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
323 $                         where first preconditioner is built from alpha I + S and second from
324 $                         alpha I + R
325 
326    Level: beginner
327 
328 .seealso: PCCompositeSetType()
329 E*/
330 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
331 PETSC_EXTERN const char *const PCCompositeTypes[];
332 
333 PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
334 PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
335 PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
336 PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
337 
338 PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
339 PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
340 PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
341 
342 PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
343 PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
344 PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
345 PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
346 PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
347 PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
348 PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
349 PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
350 
351 PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
352 PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
353 PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
354 PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
355 
356 PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
357 PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
358 PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
359 PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
360 PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
361 PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
362 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool);
363 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*);
364 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC,PetscBool);
365 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC,PetscBool*);
366 PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC,PetscBool);
367 PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC,PetscBool*);
368 
369 
370 /*E
371     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
372 
373     Level: intermediate
374 
375 .seealso: PCFieldSplitSetSchurPre()
376 E*/
377 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;
378 PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
379 
380 /*E
381     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
382 
383     Level: intermediate
384 
385 .seealso: PCFieldSplitSetSchurFactType()
386 E*/
387 typedef enum {
388   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
389   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
390   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
391   PC_FIELDSPLIT_SCHUR_FACT_FULL
392 } PCFieldSplitSchurFactType;
393 PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
394 
395 PETSC_EXTERN PETSC_DEPRECATED("Use PCFieldSplitSetSchurPre") PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
396 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC,PCFieldSplitSchurPreType,Mat);
397 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC,PCFieldSplitSchurPreType*,Mat*);
398 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
399 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
400 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC,Mat *S);
401 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC,Mat *S);
402 
403 PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
404 PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
405 
406 PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
407 
408 PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);
409 
410 PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
411 PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);
412 
413 PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
414 PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);
415 
416 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
417 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
418 PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
419 
420 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
421 PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
422 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
423 PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
424 /*E
425     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
426 
427     Level: intermediate
428 
429 .seealso: PCPARMSSetGlobal()
430 E*/
431 typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
432 PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
433 /*E
434     PCPARMSLocalType - Determines the local preconditioner method in PARMS
435 
436     Level: intermediate
437 
438 .seealso: PCPARMSSetLocal()
439 E*/
440 typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
441 PETSC_EXTERN const char *const PCPARMSLocalTypes[];
442 
443 PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC,PCPARMSGlobalType);
444 PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC,PCPARMSLocalType);
445 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC,PetscReal,PetscInt);
446 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC,PetscInt);
447 PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC,PetscBool);
448 PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC,PetscInt,PetscInt,PetscInt);
449 
450 /*E
451     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method
452 
453     Level: intermediate
454 
455 .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseProl()
456 E*/
457 typedef const char *PCGAMGType;
458 #define PCGAMGAGG         "agg"
459 #define PCGAMGGEO         "geo"
460 #define PCGAMGCLASSICAL   "classical"
461 PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
462 PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
463 PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
464 PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
465 PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
466 PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
467 PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
468 PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType);
469 PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC,PetscInt);
470 PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC,PetscBool);
471 PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool);
472 PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool);
473 PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void);
474 PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void);
475 
476 typedef const char *PCGAMGClassicalType;
477 #define PCGAMGCLASSICALDIRECT   "direct"
478 #define PCGAMGCLASSICALSTANDARD "standard"
479 PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC,PCGAMGClassicalType);
480 
481 #if defined(PETSC_HAVE_PCBDDC)
482 PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisMat(PC,Mat);
483 PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC,IS);
484 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
485 PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC,PetscInt);
486 PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace);
487 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
488 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC,IS);
489 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
490 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC,IS*);
491 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
492 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC,IS);
493 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
494 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC,IS*);
495 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
496 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC,PetscInt,IS[]);
497 PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const 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 /*E
508     PCMGType - Determines the type of multigrid method that is run.
509 
510    Level: beginner
511 
512    Values:
513 +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles()
514 .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
515                 smoothed before updating the residual. This only uses the
516                 down smoother, in the preconditioner the upper smoother is ignored
517 .  PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
518             that is starts on the coarsest grid, performs a cycle, interpolates
519             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
520             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
521             calls the V-cycle only on the coarser level and has a post-smoother instead.
522 -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
523                from a finer
524 
525 .seealso: PCMGSetType()
526 
527 E*/
528 typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
529 PETSC_EXTERN const char *const PCMGTypes[];
530 #define PC_MG_CASCADE PC_MG_KASKADE;
531 
532 /*E
533     PCMGCycleType - Use V-cycle or W-cycle
534 
535    Level: beginner
536 
537    Values:
538 +  PC_MG_V_CYCLE
539 -  PC_MG_W_CYCLE
540 
541 .seealso: PCMGSetCycleType()
542 
543 E*/
544 typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
545 PETSC_EXTERN const char *const PCMGCycleTypes[];
546 
547 PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType);
548 PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*);
549 PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*);
550 
551 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt);
552 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt);
553 PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType);
554 PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType);
555 PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt);
556 PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt);
557 PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool);
558 PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*);
559 
560 
561 PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec);
562 PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec);
563 PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec);
564 
565 PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat);
566 PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*);
567 PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat);
568 PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*);
569 PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec);
570 PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*);
571 PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
572 PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat,Vec,Vec,Vec);
573 
574 /*E
575     PCExoticType - Face based or wirebasket based coarse grid space
576 
577    Level: beginner
578 
579 .seealso: PCExoticSetType(), PCEXOTIC
580 E*/
581 typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
582 PETSC_EXTERN const char *const PCExoticTypes[];
583 PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);
584 
585 #endif /* __PETSCPC_H */
586