xref: /petsc/include/petscpc.h (revision 21e3ffae2f3b73c0bd738cf6d0a809700fc04bb0)
1 /*
2       Preconditioner module.
3 */
4 #ifndef PETSCPC_H
5 #define PETSCPC_H
6 
7 #include <petscmat.h>
8 #include <petscdmtypes.h>
9 #include <petscpctypes.h>
10 
11 /* SUBMANSEC = PC */
12 
13 PETSC_EXTERN PetscErrorCode PCInitializePackage(void);
14 
15 /*
16     PCList contains the list of preconditioners currently registered
17    These are added with PCRegister()
18 */
19 PETSC_EXTERN PetscFunctionList PCList;
20 
21 /* Logging support */
22 PETSC_EXTERN PetscClassId PC_CLASSID;
23 
24 /* Arrays of names for options in implementation PCs */
25 PETSC_EXTERN const char *const *const PCSides;
26 PETSC_EXTERN const char *const        PCJacobiTypes[];
27 PETSC_EXTERN const char *const        PCASMTypes[];
28 PETSC_EXTERN const char *const        PCGASMTypes[];
29 PETSC_EXTERN const char *const        PCCompositeTypes[];
30 PETSC_EXTERN const char *const        PCFieldSplitSchurPreTypes[];
31 PETSC_EXTERN const char *const        PCFieldSplitSchurFactTypes[];
32 PETSC_EXTERN const char *const        PCPARMSGlobalTypes[];
33 PETSC_EXTERN const char *const        PCPARMSLocalTypes[];
34 PETSC_EXTERN const char *const        PCMGTypes[];
35 PETSC_EXTERN const char *const        PCMGCycleTypes[];
36 PETSC_EXTERN const char *const        PCMGGalerkinTypes[];
37 PETSC_EXTERN const char *const        PCMGCoarseSpaceTypes[];
38 PETSC_EXTERN const char *const        PCExoticTypes[];
39 PETSC_EXTERN const char *const        PCPatchConstructTypes[];
40 PETSC_EXTERN const char *const        PCDeflationTypes[];
41 PETSC_EXTERN const char *const *const PCFailedReasons;
42 
43 PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm, PC *);
44 PETSC_EXTERN PetscErrorCode PCSetType(PC, PCType);
45 PETSC_EXTERN PetscErrorCode PCGetType(PC, PCType *);
46 PETSC_EXTERN PetscErrorCode PCSetUp(PC);
47 
48 PETSC_EXTERN PetscErrorCode PCSetFailedReason(PC, PCFailedReason);
49 PETSC_EXTERN PetscErrorCode PCGetFailedReason(PC, PCFailedReason *);
50 PETSC_DEPRECATED_FUNCTION("Use PCGetFailedReason() (since version 3.11)") static inline PetscErrorCode PCGetSetUpFailedReason(PC pc, PCFailedReason *reason)
51 {
52   return PCGetFailedReason(pc, reason);
53 }
54 PETSC_EXTERN PetscErrorCode PCGetFailedReasonRank(PC, PCFailedReason *);
55 
56 PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
57 PETSC_EXTERN PetscErrorCode PCApply(PC, Vec, Vec);
58 PETSC_EXTERN PetscErrorCode PCMatApply(PC, Mat, Mat);
59 PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC, Vec, Vec);
60 PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC, Vec, Vec);
61 PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC, PCSide, Vec, Vec, Vec);
62 PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC, Vec, Vec);
63 PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC, PetscBool *);
64 PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC, PCSide, Vec, Vec, Vec);
65 PETSC_EXTERN PetscErrorCode PCSetReusePreconditioner(PC, PetscBool);
66 PETSC_EXTERN PetscErrorCode PCGetReusePreconditioner(PC, PetscBool *);
67 PETSC_EXTERN PetscErrorCode PCSetErrorIfFailure(PC, PetscBool);
68 
69 #define PC_FILE_CLASSID 1211222
70 
71 PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *);
72 PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC, PetscBool *);
73 PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC, PetscBool);
74 PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC, PetscBool *);
75 
76 PETSC_EXTERN PetscErrorCode PCRegister(const char[], PetscErrorCode (*)(PC));
77 
78 PETSC_EXTERN PetscErrorCode PCReset(PC);
79 PETSC_EXTERN PetscErrorCode PCDestroy(PC *);
80 PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);
81 
82 PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC, Mat *);
83 PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC, PetscErrorCode (*)(PC, PetscInt, const IS[], const IS[], Mat[], void *), void *);
84 PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC, PetscInt, const IS[], const IS[], Mat[], void *);
85 
86 PETSC_EXTERN PetscErrorCode PCSetOperators(PC, Mat, Mat);
87 PETSC_EXTERN PetscErrorCode PCGetOperators(PC, Mat *, Mat *);
88 PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC, PetscBool *, PetscBool *);
89 
90 PETSC_EXTERN PetscErrorCode PCView(PC, PetscViewer);
91 PETSC_EXTERN PetscErrorCode PCLoad(PC, PetscViewer);
92 PETSC_EXTERN PetscErrorCode PCViewFromOptions(PC, PetscObject, const char[]);
93 
94 PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC, const char[]);
95 PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC, const char[]);
96 PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC, const char *[]);
97 
98 PETSC_EXTERN PetscErrorCode PCComputeOperator(PC, MatType, Mat *);
99 PETSC_DEPRECATED_FUNCTION("Use PCComputeOperator() (since version 3.12)") static inline PetscErrorCode PCComputeExplicitOperator(PC A, Mat *B)
100 {
101   return PCComputeOperator(A, NULL, B);
102 }
103 
104 /*
105       These are used to provide extra scaling of preconditioned
106    operator for time-stepping schemes like in SUNDIALS
107 */
108 PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC, PetscBool *);
109 PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC, Vec, Vec);
110 PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC, Vec, Vec);
111 PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC, Vec);
112 
113 PETSC_EXTERN PetscErrorCode PCSetDM(PC, DM);
114 PETSC_EXTERN PetscErrorCode PCGetDM(PC, DM *);
115 
116 PETSC_EXTERN PetscErrorCode PCGetInterpolations(PC, PetscInt *, Mat *[]);
117 PETSC_EXTERN PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *, Mat *[]);
118 
119 PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC, PetscInt, PetscInt, PetscReal *);
120 
121 PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC, void *);
122 PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC, void *);
123 
124 /* ------------- options specific to particular preconditioners --------- */
125 
126 PETSC_EXTERN PetscErrorCode PCJacobiSetType(PC, PCJacobiType);
127 PETSC_EXTERN PetscErrorCode PCJacobiGetType(PC, PCJacobiType *);
128 PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC, PetscBool);
129 PETSC_EXTERN PetscErrorCode PCJacobiGetUseAbs(PC, PetscBool *);
130 PETSC_EXTERN PetscErrorCode PCJacobiSetFixDiagonal(PC, PetscBool);
131 PETSC_EXTERN PetscErrorCode PCJacobiGetFixDiagonal(PC, PetscBool *);
132 PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC, MatSORType);
133 PETSC_EXTERN PetscErrorCode PCSORGetSymmetric(PC, MatSORType *);
134 PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC, PetscReal);
135 PETSC_EXTERN PetscErrorCode PCSORGetOmega(PC, PetscReal *);
136 PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC, PetscInt, PetscInt);
137 PETSC_EXTERN PetscErrorCode PCSORGetIterations(PC, PetscInt *, PetscInt *);
138 
139 PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC, PetscReal);
140 PETSC_EXTERN PetscErrorCode PCEisenstatGetOmega(PC, PetscReal *);
141 PETSC_EXTERN PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC, PetscBool);
142 PETSC_EXTERN PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC, PetscBool *);
143 
144 PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC, PetscInt, const PetscInt[]);
145 PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC, PetscInt *, const PetscInt *[]);
146 PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC, PetscInt, const PetscInt[]);
147 PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC, PetscInt *, const PetscInt *[]);
148 
149 PETSC_EXTERN PetscErrorCode PCShellSetApply(PC, PetscErrorCode (*)(PC, Vec, Vec));
150 PETSC_EXTERN PetscErrorCode PCShellSetMatApply(PC, PetscErrorCode (*)(PC, Mat, Mat));
151 PETSC_EXTERN PetscErrorCode PCShellSetApplySymmetricLeft(PC, PetscErrorCode (*)(PC, Vec, Vec));
152 PETSC_EXTERN PetscErrorCode PCShellSetApplySymmetricRight(PC, PetscErrorCode (*)(PC, Vec, Vec));
153 PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC, PetscErrorCode (*)(PC, PCSide, Vec, Vec, Vec));
154 PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC, PetscErrorCode (*)(PC, Vec, Vec));
155 PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC, PetscErrorCode (*)(PC));
156 PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC, PetscErrorCode (*)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *));
157 PETSC_EXTERN PetscErrorCode PCShellSetView(PC, PetscErrorCode (*)(PC, PetscViewer));
158 PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC, PetscErrorCode (*)(PC));
159 PETSC_EXTERN PetscErrorCode PCShellSetContext(PC, void *);
160 PETSC_EXTERN PetscErrorCode PCShellGetContext(PC, void *);
161 PETSC_EXTERN PetscErrorCode PCShellSetName(PC, const char[]);
162 PETSC_EXTERN PetscErrorCode PCShellGetName(PC, const char *[]);
163 
164 PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC, PetscReal);
165 
166 PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC, MatFactorShiftType);
167 PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC, PetscReal);
168 
169 PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverType(PC, MatSolverType);
170 PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverType(PC, MatSolverType *);
171 PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverType(PC);
172 PETSC_DEPRECATED_FUNCTION("Use PCFactorSetMatSolverType() (since version 3.9)") static inline PetscErrorCode PCFactorSetMatSolverPackage(PC pc, MatSolverType stype)
173 {
174   return PCFactorSetMatSolverType(pc, stype);
175 }
176 PETSC_DEPRECATED_FUNCTION("Use PCFactorGetMatSolverType() (since version 3.9)") static inline PetscErrorCode PCFactorGetMatSolverPackage(PC pc, MatSolverType *stype)
177 {
178   return PCFactorGetMatSolverType(pc, stype);
179 }
180 PETSC_DEPRECATED_FUNCTION("Use PCFactorSetUpMatSolverType() (since version 3.9)") static inline PetscErrorCode PCFactorSetUpMatSolverPackage(PC pc)
181 {
182   return PCFactorSetUpMatSolverType(pc);
183 }
184 
185 PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC, PetscReal);
186 PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC, PetscReal);
187 PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC, PetscReal);
188 
189 PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC, MatOrderingType);
190 PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC, PetscBool);
191 PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC, PetscBool);
192 PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC, PetscBool);
193 PETSC_EXTERN PetscErrorCode PCFactorGetUseInPlace(PC, PetscBool *);
194 PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC, PetscBool);
195 PETSC_EXTERN PetscErrorCode PCFactorGetAllowDiagonalFill(PC, PetscBool *);
196 PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC, PetscBool);
197 
198 PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC, PetscInt);
199 PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC, PetscInt *);
200 PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC, PetscReal, PetscReal, PetscInt);
201 PETSC_EXTERN PetscErrorCode PCFactorGetZeroPivot(PC, PetscReal *);
202 PETSC_EXTERN PetscErrorCode PCFactorGetShiftAmount(PC, PetscReal *);
203 PETSC_EXTERN PetscErrorCode PCFactorGetShiftType(PC, MatFactorShiftType *);
204 
205 PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC, PetscInt, IS[], IS[]);
206 PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC, PetscInt, IS[], IS[]);
207 PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC, PetscInt);
208 PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC, PetscBool);
209 PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC, PetscBool *);
210 PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC, PetscBool);
211 
212 PETSC_EXTERN PetscErrorCode PCASMSetType(PC, PCASMType);
213 PETSC_EXTERN PetscErrorCode PCASMGetType(PC, PCASMType *);
214 PETSC_EXTERN PetscErrorCode PCASMSetLocalType(PC, PCCompositeType);
215 PETSC_EXTERN PetscErrorCode PCASMGetLocalType(PC, PCCompositeType *);
216 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat, PetscInt, IS *[]);
217 PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt, IS[], IS[]);
218 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, IS **, IS **);
219 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC, PetscInt *, IS *[], IS *[]);
220 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC, PetscInt *, Mat *[]);
221 PETSC_EXTERN PetscErrorCode PCASMGetSubMatType(PC, MatType *);
222 PETSC_EXTERN PetscErrorCode PCASMSetSubMatType(PC, MatType);
223 
224 PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC, PetscInt);
225 PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC, PetscInt, IS[], IS[]);
226 PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC, PetscInt);
227 PETSC_EXTERN PetscErrorCode PCGASMSetUseDMSubdomains(PC, PetscBool);
228 PETSC_EXTERN PetscErrorCode PCGASMGetUseDMSubdomains(PC, PetscBool *);
229 PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC, PetscBool);
230 
231 PETSC_EXTERN PetscErrorCode PCGASMSetType(PC, PCGASMType);
232 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains(Mat, PetscInt, PetscInt *, IS *[]);
233 PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt, IS *[], IS *[]);
234 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, IS **, IS **);
235 PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC, PetscInt *, IS *[], IS *[]);
236 PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC, PetscInt *, Mat *[]);
237 
238 PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC, PCCompositeType);
239 PETSC_EXTERN PetscErrorCode PCCompositeGetType(PC, PCCompositeType *);
240 PETSC_EXTERN PetscErrorCode PCCompositeAddPCType(PC, PCType);
241 PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC, PC);
242 PETSC_EXTERN PetscErrorCode PCCompositeGetNumberPC(PC, PetscInt *);
243 PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC, PetscInt, PC *);
244 PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC, PetscScalar);
245 
246 PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC, PetscInt);
247 PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC, VecScatter, VecScatter);
248 PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC, Mat *, Mat *);
249 
250 PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC, PetscReal);
251 PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC, PetscInt);
252 PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC, PetscInt);
253 PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC, PetscInt);
254 PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC, PetscInt);
255 PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC, PetscInt);
256 PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC, PetscInt);
257 PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC, PetscInt);
258 
259 PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC, const char[]);
260 PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC, const char *[]);
261 PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteGradient(PC, Mat);
262 PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteCurl(PC, Mat);
263 PETSC_EXTERN PetscErrorCode PCHYPRESetInterpolations(PC, PetscInt, Mat, Mat[], Mat, Mat[]);
264 PETSC_EXTERN PetscErrorCode PCHYPRESetEdgeConstantVectors(PC, Vec, Vec, Vec);
265 PETSC_EXTERN PetscErrorCode PCHYPREAMSSetInteriorNodes(PC, Vec);
266 PETSC_EXTERN PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC, Mat);
267 PETSC_EXTERN PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC, Mat);
268 
269 PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC, const char[], PetscInt, const PetscInt *, const PetscInt *);
270 PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC, PCCompositeType);
271 PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC, PCCompositeType *);
272 PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC, PetscInt);
273 PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC, const char[], IS);
274 PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC, const char[], IS *);
275 PETSC_EXTERN PetscErrorCode PCFieldSplitGetISByIndex(PC, PetscInt, IS *);
276 PETSC_EXTERN PetscErrorCode PCFieldSplitRestrictIS(PC, IS);
277 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC, PetscBool);
278 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC, PetscBool *);
279 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC, PetscBool);
280 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC, PetscBool *);
281 PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC, PetscBool);
282 PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC, PetscBool *);
283 
284 PETSC_EXTERN                PETSC_DEPRECATED_FUNCTION("Use PCFieldSplitSetSchurPre() (since version 3.5)") PetscErrorCode PCFieldSplitSchurPrecondition(PC, PCFieldSplitSchurPreType, Mat);
285 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC, PCFieldSplitSchurPreType, Mat);
286 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC, PCFieldSplitSchurPreType *, Mat *);
287 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC, PCFieldSplitSchurFactType);
288 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurScale(PC, PetscScalar);
289 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC, Mat *, Mat *, Mat *, Mat *);
290 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC, Mat *S);
291 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC, Mat *S);
292 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDetectSaddlePoint(PC, PetscBool *);
293 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDetectSaddlePoint(PC, PetscBool);
294 PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBTol(PC, PetscReal);
295 PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBNu(PC, PetscReal);
296 PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBMaxit(PC, PetscInt);
297 PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBDelay(PC, PetscInt);
298 
299 PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC, Mat);
300 PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC, Mat);
301 PETSC_EXTERN PetscErrorCode PCGalerkinSetComputeSubmatrix(PC, PetscErrorCode (*)(PC, Mat, Mat, Mat *, void *), void *);
302 
303 PETSC_EXTERN PetscErrorCode PCPythonSetType(PC, const char[]);
304 PETSC_EXTERN PetscErrorCode PCPythonGetType(PC, const char *[]);
305 
306 PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC, PCPARMSGlobalType);
307 PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC, PCPARMSLocalType);
308 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC, PetscReal, PetscInt);
309 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC, PetscInt);
310 PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC, PetscBool);
311 PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC, PetscInt, PetscInt, PetscInt);
312 
313 PETSC_EXTERN PetscErrorCode PCGAMGSetType(PC, PCGAMGType);
314 PETSC_EXTERN PetscErrorCode PCGAMGGetType(PC, PCGAMGType *);
315 PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC, PetscInt);
316 
317 PETSC_EXTERN PetscErrorCode PCGAMGSetRepartition(PC, PetscBool);
318 PETSC_EXTERN PetscErrorCode PCGAMGSetUseSAEstEig(PC, PetscBool);
319 PETSC_EXTERN PetscErrorCode PCGAMGSetEigenvalues(PC, PetscReal, PetscReal);
320 PETSC_EXTERN PetscErrorCode PCGAMGASMSetUseAggs(PC, PetscBool);
321 PETSC_EXTERN PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC, PetscBool);
322 PETSC_EXTERN PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC, PetscBool);
323 PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC, PCGAMGLayoutType);
324 PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC, PetscReal[], PetscInt);
325 PETSC_EXTERN PetscErrorCode PCGAMGSetRankReductionFactors(PC, PetscInt[], PetscInt);
326 PETSC_EXTERN PetscErrorCode PCGAMGSetThresholdScale(PC, PetscReal);
327 PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC, PetscInt);
328 PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC, PetscInt);
329 PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC, PetscInt);
330 PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC, PetscInt);
331 PETSC_EXTERN PetscErrorCode PCGAMGSetAggressiveLevels(PC, PetscInt);
332 PETSC_EXTERN PetscErrorCode PCGAMGSetReuseInterpolation(PC, PetscBool);
333 PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void);
334 PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void);
335 PETSC_EXTERN PetscErrorCode PCGAMGRegister(PCGAMGType, PetscErrorCode (*)(PC));
336 PETSC_EXTERN PetscErrorCode PCGAMGCreateGraph(PC, Mat, Mat *);
337 
338 PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC, PCGAMGClassicalType);
339 PETSC_EXTERN PetscErrorCode PCGAMGClassicalGetType(PC, PCGAMGClassicalType *);
340 
341 PETSC_EXTERN PetscErrorCode PCBDDCSetDiscreteGradient(PC, Mat, PetscInt, PetscInt, PetscBool, PetscBool);
342 PETSC_EXTERN PetscErrorCode PCBDDCSetDivergenceMat(PC, Mat, PetscBool, IS);
343 PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisMat(PC, Mat, PetscBool);
344 PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesIS(PC, IS);
345 PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC, IS);
346 PETSC_EXTERN PetscErrorCode PCBDDCGetPrimalVerticesIS(PC, IS *);
347 PETSC_EXTERN PetscErrorCode PCBDDCGetPrimalVerticesLocalIS(PC, IS *);
348 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC, PetscInt);
349 PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC, PetscInt);
350 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC, IS);
351 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC, IS);
352 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC, IS *);
353 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC, IS *);
354 PETSC_EXTERN PetscErrorCode PCBDDCSetInterfaceExtType(PC, PCBDDCInterfaceExtType);
355 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC, IS);
356 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC, IS);
357 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC, IS *);
358 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC, IS *);
359 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC, PetscInt, IS[]);
360 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC, PetscInt, IS[]);
361 PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC, PetscInt, const PetscInt[], const PetscInt[], PetscCopyMode);
362 PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC, PetscBool, const char *, Mat *, PC *);
363 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat, Vec, Vec);
364 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat, Vec, Vec);
365 PETSC_EXTERN PetscErrorCode PCBDDCFinalizePackage(void);
366 PETSC_EXTERN PetscErrorCode PCBDDCInitializePackage(void);
367 
368 PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC, PetscBool);
369 PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC, PetscScalar);
370 PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC, Vec);
371 
372 PETSC_EXTERN PetscInt       PetscMGLevelId;
373 PETSC_EXTERN PetscErrorCode PCMGSetType(PC, PCMGType);
374 PETSC_EXTERN PetscErrorCode PCMGGetType(PC, PCMGType *);
375 PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC, PetscInt, MPI_Comm *);
376 PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC, PetscInt *);
377 
378 PETSC_EXTERN PetscErrorCode PCMGSetDistinctSmoothUp(PC);
379 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmooth(PC, PetscInt);
380 PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC, PCMGCycleType);
381 PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC, PetscInt, PCMGCycleType);
382 PETSC_DEPRECATED_FUNCTION("Use PCMGSetCycleTypeOnLevel() (since version 3.5)") static inline PetscErrorCode PCMGSetCyclesOnLevel(PC pc, PetscInt l, PetscInt t)
383 {
384   return PCMGSetCycleTypeOnLevel(pc, l, (PCMGCycleType)t);
385 }
386 PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC, PetscInt);
387 PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC, PCMGGalerkinType);
388 PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC, PCMGGalerkinType *);
389 PETSC_EXTERN PetscErrorCode PCMGSetAdaptCoarseSpaceType(PC, PCMGCoarseSpaceType);
390 PETSC_EXTERN PetscErrorCode PCMGGetAdaptCoarseSpaceType(PC, PCMGCoarseSpaceType *);
391 PETSC_EXTERN PetscErrorCode PCMGSetAdaptCR(PC, PetscBool);
392 PETSC_EXTERN PetscErrorCode PCMGGetAdaptCR(PC, PetscBool *);
393 /* MATT: Remove? */
394 PETSC_EXTERN PetscErrorCode PCMGSetAdaptInterpolation(PC, PetscBool);
395 PETSC_EXTERN PetscErrorCode PCMGGetAdaptInterpolation(PC, PetscBool *);
396 
397 PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC, PetscInt, Vec);
398 PETSC_EXTERN PetscErrorCode PCMGSetX(PC, PetscInt, Vec);
399 PETSC_EXTERN PetscErrorCode PCMGSetR(PC, PetscInt, Vec);
400 
401 PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC, PetscInt, Mat);
402 PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC, PetscInt, Mat *);
403 PETSC_EXTERN PetscErrorCode PCMGSetInjection(PC, PetscInt, Mat);
404 PETSC_EXTERN PetscErrorCode PCMGGetInjection(PC, PetscInt, Mat *);
405 PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC, PetscInt, Mat);
406 PETSC_EXTERN PetscErrorCode PCMGSetOperators(PC, PetscInt, Mat, Mat);
407 PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC, PetscInt, Mat *);
408 PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC, PetscInt, Vec);
409 PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC, PetscInt, Vec *);
410 PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC, PetscInt, PetscErrorCode (*)(Mat, Vec, Vec, Vec), Mat);
411 PETSC_EXTERN PetscErrorCode PCMGSetResidualTranspose(PC, PetscInt, PetscErrorCode (*)(Mat, Vec, Vec, Vec), Mat);
412 PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat, Vec, Vec, Vec);
413 PETSC_EXTERN PetscErrorCode PCMGResidualTransposeDefault(Mat, Vec, Vec, Vec);
414 PETSC_EXTERN PetscErrorCode PCMGMatResidualDefault(Mat, Mat, Mat, Mat);
415 PETSC_EXTERN PetscErrorCode PCMGMatResidualTransposeDefault(Mat, Mat, Mat, Mat);
416 PETSC_EXTERN PetscErrorCode PCMGGalerkinSetMatProductAlgorithm(PC, const char[]);
417 PETSC_EXTERN PetscErrorCode PCMGGalerkinGetMatProductAlgorithm(PC, const char *[]);
418 PETSC_EXTERN PetscErrorCode PCMGGetGridComplexity(PC, PetscReal *, PetscReal *);
419 
420 PETSC_EXTERN PetscErrorCode PCHMGSetReuseInterpolation(PC, PetscBool);
421 PETSC_EXTERN PetscErrorCode PCHMGSetUseSubspaceCoarsening(PC, PetscBool);
422 PETSC_EXTERN PetscErrorCode PCHMGSetInnerPCType(PC, PCType);
423 PETSC_EXTERN PetscErrorCode PCHMGSetCoarseningComponent(PC, PetscInt);
424 PETSC_EXTERN PetscErrorCode PCHMGUseMatMAIJ(PC, PetscBool);
425 
426 PETSC_EXTERN PetscErrorCode PCTelescopeGetSubcommType(PC, PetscSubcommType *);
427 PETSC_EXTERN PetscErrorCode PCTelescopeSetSubcommType(PC, PetscSubcommType);
428 PETSC_EXTERN PetscErrorCode PCTelescopeGetReductionFactor(PC, PetscInt *);
429 PETSC_EXTERN PetscErrorCode PCTelescopeSetReductionFactor(PC, PetscInt);
430 PETSC_EXTERN PetscErrorCode PCTelescopeGetIgnoreDM(PC, PetscBool *);
431 PETSC_EXTERN PetscErrorCode PCTelescopeSetIgnoreDM(PC, PetscBool);
432 PETSC_EXTERN PetscErrorCode PCTelescopeGetUseCoarseDM(PC, PetscBool *);
433 PETSC_EXTERN PetscErrorCode PCTelescopeSetUseCoarseDM(PC, PetscBool);
434 PETSC_EXTERN PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC, PetscBool *);
435 PETSC_EXTERN PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC, PetscBool);
436 PETSC_EXTERN PetscErrorCode PCTelescopeGetDM(PC, DM *);
437 
438 PETSC_EXTERN PetscErrorCode PCPatchSetSaveOperators(PC, PetscBool);
439 PETSC_EXTERN PetscErrorCode PCPatchGetSaveOperators(PC, PetscBool *);
440 PETSC_EXTERN PetscErrorCode PCPatchSetPrecomputeElementTensors(PC, PetscBool);
441 PETSC_EXTERN PetscErrorCode PCPatchGetPrecomputeElementTensors(PC, PetscBool *);
442 PETSC_EXTERN PetscErrorCode PCPatchSetPartitionOfUnity(PC, PetscBool);
443 PETSC_EXTERN PetscErrorCode PCPatchGetPartitionOfUnity(PC, PetscBool *);
444 PETSC_EXTERN PetscErrorCode PCPatchSetSubMatType(PC, MatType);
445 PETSC_EXTERN PetscErrorCode PCPatchGetSubMatType(PC, MatType *);
446 PETSC_EXTERN PetscErrorCode PCPatchSetCellNumbering(PC, PetscSection);
447 PETSC_EXTERN PetscErrorCode PCPatchGetCellNumbering(PC, PetscSection *);
448 PETSC_EXTERN PetscErrorCode PCPatchSetConstructType(PC, PCPatchConstructType, PetscErrorCode (*)(PC, PetscInt *, IS **, IS *, void *), void *);
449 PETSC_EXTERN PetscErrorCode PCPatchGetConstructType(PC, PCPatchConstructType *, PetscErrorCode (**)(PC, PetscInt *, IS **, IS *, void *), void **);
450 PETSC_EXTERN PetscErrorCode PCPatchSetDiscretisationInfo(PC, PetscInt, DM *, PetscInt *, PetscInt *, const PetscInt **, const PetscInt *, PetscInt, const PetscInt *, PetscInt, const PetscInt *);
451 PETSC_EXTERN PetscErrorCode PCPatchSetComputeOperator(PC, PetscErrorCode (*)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
452 PETSC_EXTERN PetscErrorCode PCPatchSetComputeFunction(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx);
453 PETSC_EXTERN PetscErrorCode PCPatchSetComputeOperatorInteriorFacets(PC, PetscErrorCode (*)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
454 PETSC_EXTERN PetscErrorCode PCPatchSetComputeFunctionInteriorFacets(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx);
455 
456 PETSC_EXTERN PetscErrorCode PCLMVMSetMatLMVM(PC, Mat);
457 PETSC_EXTERN PetscErrorCode PCLMVMGetMatLMVM(PC, Mat *);
458 PETSC_EXTERN PetscErrorCode PCLMVMSetIS(PC, IS);
459 PETSC_EXTERN PetscErrorCode PCLMVMClearIS(PC);
460 
461 PETSC_EXTERN PetscErrorCode PCExoticSetType(PC, PCExoticType);
462 
463 PETSC_EXTERN PetscErrorCode PCDeflationSetInitOnly(PC, PetscBool);
464 PETSC_EXTERN PetscErrorCode PCDeflationSetLevels(PC, PetscInt);
465 PETSC_EXTERN PetscErrorCode PCDeflationSetReductionFactor(PC, PetscInt);
466 PETSC_EXTERN PetscErrorCode PCDeflationSetCorrectionFactor(PC, PetscScalar);
467 PETSC_EXTERN PetscErrorCode PCDeflationSetSpaceToCompute(PC, PCDeflationSpaceType, PetscInt);
468 PETSC_EXTERN PetscErrorCode PCDeflationSetSpace(PC, Mat, PetscBool);
469 PETSC_EXTERN PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC, Mat);
470 PETSC_EXTERN PetscErrorCode PCDeflationSetCoarseMat(PC, Mat);
471 PETSC_EXTERN PetscErrorCode PCDeflationGetPC(PC, PC *);
472 
473 PETSC_EXTERN PetscErrorCode PCHPDDMSetAuxiliaryMat(PC, IS, Mat, PetscErrorCode (*)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void *), void *);
474 PETSC_EXTERN PetscErrorCode PCHPDDMSetRHSMat(PC, Mat);
475 PETSC_EXTERN PetscErrorCode PCHPDDMHasNeumannMat(PC, PetscBool);
476 PETSC_EXTERN PetscErrorCode PCHPDDMSetCoarseCorrectionType(PC, PCHPDDMCoarseCorrectionType);
477 PETSC_EXTERN PetscErrorCode PCHPDDMGetCoarseCorrectionType(PC, PCHPDDMCoarseCorrectionType *);
478 PETSC_EXTERN PetscErrorCode PCHPDDMSetSTShareSubKSP(PC, PetscBool);
479 PETSC_EXTERN PetscErrorCode PCHPDDMGetSTShareSubKSP(PC, PetscBool *);
480 PETSC_EXTERN PetscErrorCode PCHPDDMSetDeflationMat(PC, IS, Mat);
481 PETSC_EXTERN PetscErrorCode PCHPDDMFinalizePackage(void);
482 PETSC_EXTERN PetscErrorCode PCHPDDMInitializePackage(void);
483 
484 PETSC_EXTERN PetscErrorCode PCAmgXGetResources(PC, void *);
485 
486 #endif /* PETSCPC_H */
487