Actual source code: petscpc.h

  1: /*
  2:       Preconditioner module.
  3: */
  4: #pragma once

  6: #include <petscmat.h>
  7: #include <petscdmtypes.h>
  8: #include <petscpctypes.h>

 10: /* MANSEC = KSP */
 11: /* SUBMANSEC = PC */

 13: PETSC_EXTERN PetscErrorCode PCInitializePackage(void);
 14: PETSC_EXTERN PetscErrorCode PCFinalizePackage(void);

 16: /*
 17:     PCList contains the list of preconditioners currently registered
 18:    These are added with PCRegister()
 19: */
 20: PETSC_EXTERN PetscFunctionList PCList;

 22: /* Logging support */
 23: PETSC_EXTERN PetscClassId PC_CLASSID;

 25: /* Arrays of names for options in implementation PCs */
 26: PETSC_EXTERN const char *const *const PCSides;
 27: PETSC_EXTERN const char *const        PCJacobiTypes[];
 28: PETSC_EXTERN const char *const        PCASMTypes[];
 29: PETSC_EXTERN const char *const        PCGASMTypes[];
 30: PETSC_EXTERN const char *const        PCCompositeTypes[];
 31: PETSC_EXTERN const char *const        PCFieldSplitSchurPreTypes[];
 32: PETSC_EXTERN const char *const        PCFieldSplitSchurFactTypes[];
 33: PETSC_EXTERN const char *const        PCPARMSGlobalTypes[];
 34: PETSC_EXTERN const char *const        PCPARMSLocalTypes[];
 35: PETSC_EXTERN const char *const        PCMGTypes[];
 36: PETSC_EXTERN const char *const        PCMGCycleTypes[];
 37: PETSC_EXTERN const char *const        PCMGGalerkinTypes[];
 38: PETSC_EXTERN const char *const        PCMGCoarseSpaceTypes[];
 39: PETSC_EXTERN const char *const        PCExoticTypes[];
 40: PETSC_EXTERN const char *const        PCPatchConstructTypes[];
 41: PETSC_EXTERN const char *const        PCDeflationTypes[];
 42: PETSC_EXTERN const char *const *const PCFailedReasons;

 44: PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm, PC *);
 45: PETSC_EXTERN PetscErrorCode PCSetType(PC, PCType);
 46: PETSC_EXTERN PetscErrorCode PCGetType(PC, PCType *);
 47: PETSC_EXTERN PetscErrorCode PCSetUp(PC);

 49: PETSC_EXTERN PetscErrorCode PCSetKSPNestLevel(PC, PetscInt);
 50: PETSC_EXTERN PetscErrorCode PCGetKSPNestLevel(PC, PetscInt *);

 52: PETSC_EXTERN PetscErrorCode PCSetFailedReason(PC, PCFailedReason);
 53: PETSC_EXTERN PetscErrorCode PCGetFailedReason(PC, PCFailedReason *);
 54: PETSC_DEPRECATED_FUNCTION(3, 11, 0, "PCGetFailedReason()", ) static inline PetscErrorCode PCGetSetUpFailedReason(PC pc, PCFailedReason *reason)
 55: {
 56:   return PCGetFailedReason(pc, reason);
 57: }
 58: PETSC_EXTERN PetscErrorCode PCReduceFailedReason(PC);

 60: PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
 61: PETSC_EXTERN PetscErrorCode PCApply(PC, Vec, Vec);
 62: PETSC_EXTERN PetscErrorCode PCMatApply(PC, Mat, Mat);
 63: PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC, Vec, Vec);
 64: PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC, Vec, Vec);
 65: PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC, PCSide, Vec, Vec, Vec);
 66: PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC, Vec, Vec);
 67: PETSC_EXTERN PetscErrorCode PCMatApplyTranspose(PC, Mat, Mat);
 68: PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC, PetscBool *);
 69: PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC, PCSide, Vec, Vec, Vec);
 70: PETSC_EXTERN PetscErrorCode PCSetReusePreconditioner(PC, PetscBool);
 71: PETSC_EXTERN PetscErrorCode PCGetReusePreconditioner(PC, PetscBool *);
 72: PETSC_EXTERN PetscErrorCode PCSetErrorIfFailure(PC, PetscBool);

 74: #define PC_FILE_CLASSID 1211222

 76: PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *);
 77: PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC, PetscBool *);
 78: PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC, PetscBool);
 79: PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC, PetscBool *);

 81: PETSC_EXTERN PetscErrorCode PCRegister(const char[], PetscErrorCode (*)(PC));

 83: PETSC_EXTERN PetscErrorCode PCReset(PC);
 84: PETSC_EXTERN PetscErrorCode PCDestroy(PC *);
 85: PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);

 87: PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC, Mat *);

 89: /*S
 90:   PCModifySubMatricesFn - A prototype of a function used to modify submatrices generated with `PCASM`, `PCBJACOBI`, etc.

 92:   Calling Sequence:
 93: + pc     - the `PC` preconditioner context
 94: . nsub   - number of index sets
 95: . row    - an array of index sets that contain the global row numbers
 96:          that comprise each local submatrix
 97: . col    - an array of index sets that contain the global column numbers
 98:          that comprise each local submatrix
 99: . submat - array of local submatrices
100: - ctx    - optional user-defined context for private data for the user-defined func routine (may be `NULL`), provided with `PCSetModifySubMatrices()`

102:   Level: beginner

104: .seealso: [](ch_ksp), `PC`, `PCSetModifySubMatrices()`, `PCModifySubMatrices()`, `PCASM`, `PCBJACOBI`
105: S*/
106: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PCModifySubMatricesFn(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx);

108: PETSC_EXTERN PetscErrorCode        PCSetModifySubMatrices(PC, PCModifySubMatricesFn *, void *);
109: PETSC_EXTERN PCModifySubMatricesFn PCModifySubMatrices;

111: PETSC_EXTERN PetscErrorCode PCSetOperators(PC, Mat, Mat);
112: PETSC_EXTERN PetscErrorCode PCGetOperators(PC, Mat *, Mat *);
113: PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC, PetscBool *, PetscBool *);

115: PETSC_EXTERN PetscErrorCode PCView(PC, PetscViewer);
116: PETSC_EXTERN PetscErrorCode PCLoad(PC, PetscViewer);
117: PETSC_EXTERN PetscErrorCode PCViewFromOptions(PC, PetscObject, const char[]);

119: PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC, const char[]);
120: PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC, const char[]);
121: PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC, const char *[]);

123: PETSC_EXTERN PetscErrorCode PCComputeOperator(PC, MatType, Mat *);
124: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "PCComputeOperator()", ) static inline PetscErrorCode PCComputeExplicitOperator(PC A, Mat *B)
125: {
126:   return PCComputeOperator(A, PETSC_NULLPTR, B);
127: }

129: PETSC_EXTERN PetscErrorCode PCSetPostSetUp(PC, PetscErrorCode (*)(PC));

131: /*
132:       These are used to provide extra scaling of preconditioned
133:    operator for time-stepping schemes like in SUNDIALS
134: */
135: PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC, PetscBool *);
136: PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC, Vec, Vec);
137: PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC, Vec, Vec);
138: PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC, Vec);

140: PETSC_EXTERN PetscErrorCode PCSetDM(PC, DM);
141: PETSC_EXTERN PetscErrorCode PCGetDM(PC, DM *);

143: PETSC_EXTERN PetscErrorCode PCGetInterpolations(PC, PetscInt *, Mat *[]);
144: PETSC_EXTERN PetscErrorCode PCGetCoarseOperators(PC, PetscInt *, Mat *[]);

146: PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC, PetscInt, PetscInt, PetscReal[]);

148: PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC, void *);
149: PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC, void *);

151: /* ------------- options specific to particular preconditioners --------- */

153: PETSC_EXTERN PetscErrorCode PCJacobiSetType(PC, PCJacobiType);
154: PETSC_EXTERN PetscErrorCode PCJacobiGetType(PC, PCJacobiType *);
155: PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC, PetscBool);
156: PETSC_EXTERN PetscErrorCode PCJacobiGetUseAbs(PC, PetscBool *);
157: PETSC_EXTERN PetscErrorCode PCJacobiSetFixDiagonal(PC, PetscBool);
158: PETSC_EXTERN PetscErrorCode PCJacobiGetFixDiagonal(PC, PetscBool *);
159: PETSC_EXTERN PetscErrorCode PCJacobiGetDiagonal(PC, Vec, Vec);
160: PETSC_EXTERN PetscErrorCode PCJacobiSetRowl1Scale(PC, PetscReal);
161: PETSC_EXTERN PetscErrorCode PCJacobiGetRowl1Scale(PC, PetscReal *);
162: PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC, MatSORType);
163: PETSC_EXTERN PetscErrorCode PCSORGetSymmetric(PC, MatSORType *);
164: PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC, PetscReal);
165: PETSC_EXTERN PetscErrorCode PCSORGetOmega(PC, PetscReal *);
166: PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC, PetscInt, PetscInt);
167: PETSC_EXTERN PetscErrorCode PCSORGetIterations(PC, PetscInt *, PetscInt *);

169: PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC, PetscReal);
170: PETSC_EXTERN PetscErrorCode PCEisenstatGetOmega(PC, PetscReal *);
171: PETSC_EXTERN PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC, PetscBool);
172: PETSC_EXTERN PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC, PetscBool *);

174: PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC, PetscInt, const PetscInt[]);
175: PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC, PetscInt *, const PetscInt *[]);
176: PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC, PetscInt, const PetscInt[]);
177: PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC, PetscInt *, const PetscInt *[]);

179: PETSC_EXTERN PetscErrorCode PCShellSetApply(PC, PetscErrorCode (*)(PC, Vec, Vec));
180: PETSC_EXTERN PetscErrorCode PCShellSetMatApply(PC, PetscErrorCode (*)(PC, Mat, Mat));
181: PETSC_EXTERN PetscErrorCode PCShellSetApplySymmetricLeft(PC, PetscErrorCode (*)(PC, Vec, Vec));
182: PETSC_EXTERN PetscErrorCode PCShellSetApplySymmetricRight(PC, PetscErrorCode (*)(PC, Vec, Vec));
183: PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC, PetscErrorCode (*)(PC, PCSide, Vec, Vec, Vec));
184: PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC, PetscErrorCode (*)(PC, Vec, Vec));
185: PETSC_EXTERN PetscErrorCode PCShellSetMatApplyTranspose(PC, PetscErrorCode (*)(PC, Mat, Mat));
186: PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC, PetscErrorCode (*)(PC));
187: PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC, PetscErrorCode (*)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *));
188: PETSC_EXTERN PetscErrorCode PCShellSetView(PC, PetscErrorCode (*)(PC, PetscViewer));
189: PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC, PetscErrorCode (*)(PC));
190: PETSC_EXTERN PetscErrorCode PCShellSetContext(PC, void *);
191: PETSC_EXTERN PetscErrorCode PCShellGetContext(PC, void *);
192: PETSC_EXTERN PetscErrorCode PCShellSetName(PC, const char[]);
193: PETSC_EXTERN PetscErrorCode PCShellGetName(PC, const char *[]);

195: PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC, PetscReal);

197: PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC, MatFactorShiftType);
198: PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC, PetscReal);

200: PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverType(PC, MatSolverType);
201: PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverType(PC, MatSolverType *);
202: PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverType(PC);
203: PETSC_DEPRECATED_FUNCTION(3, 9, 0, "PCFactorSetMatSolverType()", ) static inline PetscErrorCode PCFactorSetMatSolverPackage(PC pc, MatSolverType stype)
204: {
205:   return PCFactorSetMatSolverType(pc, stype);
206: }
207: PETSC_DEPRECATED_FUNCTION(3, 9, 0, "PCFactorGetMatSolverType()", ) static inline PetscErrorCode PCFactorGetMatSolverPackage(PC pc, MatSolverType *stype)
208: {
209:   return PCFactorGetMatSolverType(pc, stype);
210: }
211: PETSC_DEPRECATED_FUNCTION(3, 9, 0, "PCFactorSetUpMatSolverType()", ) static inline PetscErrorCode PCFactorSetUpMatSolverPackage(PC pc)
212: {
213:   return PCFactorSetUpMatSolverType(pc);
214: }

216: PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC, PetscReal);
217: PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC, PetscReal);
218: PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC, PetscReal);

220: PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC, MatOrderingType);
221: PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC, PetscBool);
222: PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC, PetscBool);
223: PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC, PetscBool);
224: PETSC_EXTERN PetscErrorCode PCFactorGetUseInPlace(PC, PetscBool *);
225: PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC, PetscBool);
226: PETSC_EXTERN PetscErrorCode PCFactorGetAllowDiagonalFill(PC, PetscBool *);
227: PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC, PetscBool);

229: PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC, PetscInt);
230: PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC, PetscInt *);
231: PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC, PetscReal, PetscReal, PetscInt);
232: PETSC_EXTERN PetscErrorCode PCFactorGetZeroPivot(PC, PetscReal *);
233: PETSC_EXTERN PetscErrorCode PCFactorGetShiftAmount(PC, PetscReal *);
234: PETSC_EXTERN PetscErrorCode PCFactorGetShiftType(PC, MatFactorShiftType *);

236: PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC, PetscInt, IS[], IS[]);
237: PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC, PetscInt, IS[], IS[]);
238: PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC, PetscInt);
239: PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC, PetscBool);
240: PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC, PetscBool *);
241: PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC, PetscBool);

243: PETSC_EXTERN PetscErrorCode PCASMSetType(PC, PCASMType);
244: PETSC_EXTERN PetscErrorCode PCASMGetType(PC, PCASMType *);
245: PETSC_EXTERN PetscErrorCode PCASMSetLocalType(PC, PCCompositeType);
246: PETSC_EXTERN PetscErrorCode PCASMGetLocalType(PC, PCCompositeType *);
247: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat, PetscInt, IS *[]);
248: PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt, IS *[], IS *[]);
249: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, IS *[], IS *[]);
250: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC, PetscInt *, IS *[], IS *[]);
251: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC, PetscInt *, Mat *[]);
252: PETSC_EXTERN PetscErrorCode PCASMGetSubMatType(PC, MatType *);
253: PETSC_EXTERN PetscErrorCode PCASMSetSubMatType(PC, MatType);

255: PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC, PetscInt);
256: PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC, PetscInt, IS[], IS[]);
257: PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC, PetscInt);
258: PETSC_EXTERN PetscErrorCode PCGASMSetUseDMSubdomains(PC, PetscBool);
259: PETSC_EXTERN PetscErrorCode PCGASMGetUseDMSubdomains(PC, PetscBool *);
260: PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC, PetscBool);

262: PETSC_EXTERN PetscErrorCode PCGASMSetType(PC, PCGASMType);
263: PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains(Mat, PetscInt, PetscInt *, IS *[]);
264: PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt, IS *[], IS *[]);
265: PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, IS *[], IS *[]);
266: PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC, PetscInt *, IS *[], IS *[]);
267: PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC, PetscInt *, Mat *[]);

269: PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC, PCCompositeType);
270: PETSC_EXTERN PetscErrorCode PCCompositeGetType(PC, PCCompositeType *);
271: PETSC_EXTERN PetscErrorCode PCCompositeAddPCType(PC, PCType);
272: PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC, PC);
273: PETSC_EXTERN PetscErrorCode PCCompositeGetNumberPC(PC, PetscInt *);
274: PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC, PetscInt, PC *);
275: PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC, PetscScalar);
276: PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlphaMat(PC, Mat);

278: PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC, PetscInt);
279: PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC, VecScatter, VecScatter);
280: PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC, Mat *, Mat *);

282: PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC, PetscReal);
283: PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC, PetscInt);
284: PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC, PetscInt);
285: PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC, PetscInt);
286: PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC, PetscInt);
287: PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC, PetscInt);
288: PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC, PetscInt);
289: PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC, PetscInt);

291: PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC, const char[]);
292: PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC, const char *[]);
293: PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteGradient(PC, Mat);
294: PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteCurl(PC, Mat);
295: PETSC_EXTERN PetscErrorCode PCHYPRESetInterpolations(PC, PetscInt, Mat, Mat[], Mat, Mat[]);
296: PETSC_EXTERN PetscErrorCode PCHYPRESetEdgeConstantVectors(PC, Vec, Vec, Vec);
297: PETSC_EXTERN PetscErrorCode PCHYPREAMSSetInteriorNodes(PC, Vec);
298: PETSC_EXTERN PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC, Mat);
299: PETSC_EXTERN PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC, Mat);
300: PETSC_EXTERN PetscErrorCode PCHYPREGetCFMarkers(PC pc, PetscInt *[], PetscBT *[]);

302: PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC, const char[], PetscInt, const PetscInt *, const PetscInt *);
303: PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC, PCCompositeType);
304: PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC, PCCompositeType *);
305: PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC, PetscInt);
306: PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC, const char[], IS);
307: PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC, const char[], IS *);
308: PETSC_EXTERN PetscErrorCode PCFieldSplitGetISByIndex(PC, PetscInt, IS *);
309: PETSC_EXTERN PetscErrorCode PCFieldSplitRestrictIS(PC, IS);
310: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC, PetscBool);
311: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC, PetscBool *);
312: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC, PetscBool);
313: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC, PetscBool *);
314: PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC, PetscBool);
315: PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC, PetscBool *);

317: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "PCFieldSplitSetSchurPre()", ) PetscErrorCode PCFieldSplitSchurPrecondition(PC, PCFieldSplitSchurPreType, Mat);
318: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC, PCFieldSplitSchurPreType, Mat);
319: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC, PCFieldSplitSchurPreType *, Mat *);
320: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC, PCFieldSplitSchurFactType);
321: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurScale(PC, PetscScalar);
322: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC, Mat *, Mat *, Mat *, Mat *);
323: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC, Mat *);
324: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC, Mat *);
325: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDetectSaddlePoint(PC, PetscBool *);
326: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDetectSaddlePoint(PC, PetscBool);
327: PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBTol(PC, PetscReal);
328: PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBNu(PC, PetscReal);
329: PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBMaxit(PC, PetscInt);
330: PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBDelay(PC, PetscInt);

332: PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC, Mat);
333: PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC, Mat);
334: PETSC_EXTERN PetscErrorCode PCGalerkinSetComputeSubmatrix(PC, PetscErrorCode (*)(PC, Mat, Mat, Mat *, void *), void *);

336: PETSC_EXTERN PetscErrorCode PCPythonSetType(PC, const char[]);
337: PETSC_EXTERN PetscErrorCode PCPythonGetType(PC, const char *[]);

339: PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC, PCPARMSGlobalType);
340: PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC, PCPARMSLocalType);
341: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC, PetscReal, PetscInt);
342: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC, PetscInt);
343: PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC, PetscBool);
344: PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC, PetscInt, PetscInt, PetscInt);

346: PETSC_EXTERN PetscErrorCode PCGAMGSetType(PC, PCGAMGType);
347: PETSC_EXTERN PetscErrorCode PCGAMGGetType(PC, PCGAMGType *);
348: PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC, PetscInt);

350: PETSC_EXTERN PetscErrorCode PCGAMGSetRepartition(PC, PetscBool);
351: PETSC_EXTERN PetscErrorCode PCGAMGSetUseSAEstEig(PC, PetscBool);
352: PETSC_EXTERN PetscErrorCode PCGAMGSetRecomputeEstEig(PC, PetscBool);
353: PETSC_EXTERN PetscErrorCode PCGAMGSetEigenvalues(PC, PetscReal, PetscReal);
354: PETSC_EXTERN PetscErrorCode PCGAMGASMSetUseAggs(PC, PetscBool);
355: PETSC_EXTERN PetscErrorCode PCGAMGSetParallelCoarseGridSolve(PC, PetscBool);
356: PETSC_EXTERN PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC, PetscBool);
357: PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC, PCGAMGLayoutType);
358: PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC, PetscReal[], PetscInt);
359: PETSC_EXTERN PetscErrorCode PCGAMGSetRankReductionFactors(PC, PetscInt[], PetscInt);
360: PETSC_EXTERN PetscErrorCode PCGAMGSetThresholdScale(PC, PetscReal);
361: PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC, PetscInt);
362: PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC, PetscInt);
363: PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC, PetscInt);
364: PETSC_EXTERN PetscErrorCode PCGAMGSetAggressiveLevels(PC, PetscInt);
365: PETSC_EXTERN PetscErrorCode PCGAMGSetReuseInterpolation(PC, PetscBool);
366: PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void);
367: PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void);
368: PETSC_EXTERN PetscErrorCode PCGAMGRegister(PCGAMGType, PetscErrorCode (*)(PC));
369: PETSC_EXTERN PetscErrorCode PCGAMGCreateGraph(PC, Mat, Mat *);
370: PETSC_EXTERN PetscErrorCode PCGAMGSetAggressiveSquareGraph(PC, PetscBool);
371: PETSC_EXTERN PetscErrorCode PCGAMGMISkSetMinDegreeOrdering(PC, PetscBool);
372: PETSC_EXTERN PetscErrorCode PCGAMGMISkSetAggressive(PC, PetscInt);
373: PETSC_EXTERN PetscErrorCode PCGAMGASMSetHEM(PC, PetscInt);
374: PETSC_EXTERN PetscErrorCode PCGAMGSetLowMemoryFilter(PC, PetscBool);
375: PETSC_EXTERN PetscErrorCode PCGAMGSetGraphSymmetrize(PC, PetscBool);
376: PETSC_EXTERN PetscErrorCode PCGAMGSetInjectionIndex(PC, PetscInt, PetscInt[]);

378: PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC, PCGAMGClassicalType);
379: PETSC_EXTERN PetscErrorCode PCGAMGClassicalGetType(PC, PCGAMGClassicalType *);

381: PETSC_EXTERN PetscErrorCode PCBDDCSetDiscreteGradient(PC, Mat, PetscInt, PetscInt, PetscBool, PetscBool);
382: PETSC_EXTERN PetscErrorCode PCBDDCSetDivergenceMat(PC, Mat, PetscBool, IS);
383: PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisMat(PC, Mat, PetscBool);
384: PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesIS(PC, IS);
385: PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC, IS);
386: PETSC_EXTERN PetscErrorCode PCBDDCGetPrimalVerticesIS(PC, IS *);
387: PETSC_EXTERN PetscErrorCode PCBDDCGetPrimalVerticesLocalIS(PC, IS *);
388: PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC, PetscInt);
389: PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC, PetscInt);
390: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC, IS);
391: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC, IS);
392: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC, IS *);
393: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC, IS *);
394: PETSC_EXTERN PetscErrorCode PCBDDCSetInterfaceExtType(PC, PCBDDCInterfaceExtType);
395: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC, IS);
396: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC, IS);
397: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC, IS *);
398: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC, IS *);
399: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC, PetscInt, IS[]);
400: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC, PetscInt, IS[]);
401: PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC, PetscInt, const PetscInt[], const PetscInt[], PetscCopyMode);
402: PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC, PetscBool, const char *, Mat *, PC *);
403: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat, Vec, Vec);
404: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat, Vec, Vec);
405: PETSC_EXTERN PetscErrorCode PCBDDCFinalizePackage(void);
406: PETSC_EXTERN PetscErrorCode PCBDDCInitializePackage(void);

408: PETSC_EXTERN PetscErrorCode PCISInitialize(PC);
409: PETSC_EXTERN PetscErrorCode PCISSetUp(PC, PetscBool, PetscBool);
410: PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC, PetscBool);
411: PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC, PetscScalar);
412: PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC, Vec);
413: PETSC_EXTERN PetscErrorCode PCISScatterArrayNToVecB(PC, PetscScalar *, Vec, InsertMode, ScatterMode);
414: PETSC_EXTERN PetscErrorCode PCISApplySchur(PC, Vec, Vec, Vec, Vec, Vec);
415: PETSC_EXTERN PetscErrorCode PCISApplyInvSchur(PC, Vec, Vec, Vec, Vec);
416: PETSC_EXTERN PetscErrorCode PCISReset(PC);

418: PETSC_EXTERN PetscInt       PetscMGLevelId;
419: PETSC_EXTERN PetscErrorCode PCMGSetType(PC, PCMGType);
420: PETSC_EXTERN PetscErrorCode PCMGGetType(PC, PCMGType *);
421: PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC, PetscInt, MPI_Comm *);
422: PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC, PetscInt *);

424: PETSC_EXTERN PetscErrorCode PCMGSetDistinctSmoothUp(PC);
425: PETSC_EXTERN PetscErrorCode PCMGSetNumberSmooth(PC, PetscInt);
426: PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC, PCMGCycleType);
427: PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC, PetscInt, PCMGCycleType);
428: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "PCMGSetCycleTypeOnLevel()", ) static inline PetscErrorCode PCMGSetCyclesOnLevel(PC pc, PetscInt l, PetscInt t)
429: {
430:   return PCMGSetCycleTypeOnLevel(pc, l, (PCMGCycleType)t);
431: }
432: PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC, PetscInt);
433: PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC, PCMGGalerkinType);
434: PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC, PCMGGalerkinType *);
435: PETSC_EXTERN PetscErrorCode PCMGSetAdaptCoarseSpaceType(PC, PCMGCoarseSpaceType);
436: PETSC_EXTERN PetscErrorCode PCMGGetAdaptCoarseSpaceType(PC, PCMGCoarseSpaceType *);
437: PETSC_EXTERN PetscErrorCode PCMGSetAdaptCR(PC, PetscBool);
438: PETSC_EXTERN PetscErrorCode PCMGGetAdaptCR(PC, PetscBool *);
439: /* MATT: Remove? */
440: PETSC_EXTERN PetscErrorCode PCMGSetAdaptInterpolation(PC, PetscBool);
441: PETSC_EXTERN PetscErrorCode PCMGGetAdaptInterpolation(PC, PetscBool *);

443: PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC, PetscInt, Vec);
444: PETSC_EXTERN PetscErrorCode PCMGSetX(PC, PetscInt, Vec);
445: PETSC_EXTERN PetscErrorCode PCMGSetR(PC, PetscInt, Vec);

447: PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC, PetscInt, Mat);
448: PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC, PetscInt, Mat *);
449: PETSC_EXTERN PetscErrorCode PCMGSetInjection(PC, PetscInt, Mat);
450: PETSC_EXTERN PetscErrorCode PCMGGetInjection(PC, PetscInt, Mat *);
451: PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC, PetscInt, Mat);
452: PETSC_EXTERN PetscErrorCode PCMGSetOperators(PC, PetscInt, Mat, Mat);
453: PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC, PetscInt, Mat *);
454: PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC, PetscInt, Vec);
455: PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC, PetscInt, Vec *);
456: PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC, PetscInt, PetscErrorCode (*)(Mat, Vec, Vec, Vec), Mat);
457: PETSC_EXTERN PetscErrorCode PCMGSetResidualTranspose(PC, PetscInt, PetscErrorCode (*)(Mat, Vec, Vec, Vec), Mat);
458: PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat, Vec, Vec, Vec);
459: PETSC_EXTERN PetscErrorCode PCMGResidualTransposeDefault(Mat, Vec, Vec, Vec);
460: PETSC_EXTERN PetscErrorCode PCMGMatResidualDefault(Mat, Mat, Mat, Mat);
461: PETSC_EXTERN PetscErrorCode PCMGMatResidualTransposeDefault(Mat, Mat, Mat, Mat);
462: PETSC_EXTERN PetscErrorCode PCMGGalerkinSetMatProductAlgorithm(PC, const char[]);
463: PETSC_EXTERN PetscErrorCode PCMGGalerkinGetMatProductAlgorithm(PC, const char *[]);
464: PETSC_EXTERN PetscErrorCode PCMGGetGridComplexity(PC, PetscReal *, PetscReal *);

466: PETSC_EXTERN PetscErrorCode PCHMGSetReuseInterpolation(PC, PetscBool);
467: PETSC_EXTERN PetscErrorCode PCHMGSetUseSubspaceCoarsening(PC, PetscBool);
468: PETSC_EXTERN PetscErrorCode PCHMGSetInnerPCType(PC, PCType);
469: PETSC_EXTERN PetscErrorCode PCHMGSetCoarseningComponent(PC, PetscInt);
470: PETSC_EXTERN PetscErrorCode PCHMGUseMatMAIJ(PC, PetscBool);

472: PETSC_EXTERN PetscErrorCode PCTelescopeGetSubcommType(PC, PetscSubcommType *);
473: PETSC_EXTERN PetscErrorCode PCTelescopeSetSubcommType(PC, PetscSubcommType);
474: PETSC_EXTERN PetscErrorCode PCTelescopeGetReductionFactor(PC, PetscInt *);
475: PETSC_EXTERN PetscErrorCode PCTelescopeSetReductionFactor(PC, PetscInt);
476: PETSC_EXTERN PetscErrorCode PCTelescopeGetIgnoreDM(PC, PetscBool *);
477: PETSC_EXTERN PetscErrorCode PCTelescopeSetIgnoreDM(PC, PetscBool);
478: PETSC_EXTERN PetscErrorCode PCTelescopeGetUseCoarseDM(PC, PetscBool *);
479: PETSC_EXTERN PetscErrorCode PCTelescopeSetUseCoarseDM(PC, PetscBool);
480: PETSC_EXTERN PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC, PetscBool *);
481: PETSC_EXTERN PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC, PetscBool);
482: PETSC_EXTERN PetscErrorCode PCTelescopeGetDM(PC, DM *);

484: PETSC_EXTERN PetscErrorCode PCPatchSetSaveOperators(PC, PetscBool);
485: PETSC_EXTERN PetscErrorCode PCPatchGetSaveOperators(PC, PetscBool *);
486: PETSC_EXTERN PetscErrorCode PCPatchSetPrecomputeElementTensors(PC, PetscBool);
487: PETSC_EXTERN PetscErrorCode PCPatchGetPrecomputeElementTensors(PC, PetscBool *);
488: PETSC_EXTERN PetscErrorCode PCPatchSetPartitionOfUnity(PC, PetscBool);
489: PETSC_EXTERN PetscErrorCode PCPatchGetPartitionOfUnity(PC, PetscBool *);
490: PETSC_EXTERN PetscErrorCode PCPatchSetSubMatType(PC, MatType);
491: PETSC_EXTERN PetscErrorCode PCPatchGetSubMatType(PC, MatType *);
492: PETSC_EXTERN PetscErrorCode PCPatchSetCellNumbering(PC, PetscSection);
493: PETSC_EXTERN PetscErrorCode PCPatchGetCellNumbering(PC, PetscSection *);
494: PETSC_EXTERN PetscErrorCode PCPatchSetConstructType(PC, PCPatchConstructType, PetscErrorCode (*)(PC, PetscInt *, IS **, IS *, void *), void *);
495: PETSC_EXTERN PetscErrorCode PCPatchGetConstructType(PC, PCPatchConstructType *, PetscErrorCode (**)(PC, PetscInt *, IS **, IS *, void *), void **);
496: PETSC_EXTERN PetscErrorCode PCPatchSetDiscretisationInfo(PC, PetscInt, DM *, PetscInt *, PetscInt *, const PetscInt **, const PetscInt *, PetscInt, const PetscInt *, PetscInt, const PetscInt *);
497: PETSC_EXTERN PetscErrorCode PCPatchSetComputeOperator(PC, PetscErrorCode (*)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
498: PETSC_EXTERN PetscErrorCode PCPatchSetComputeFunction(PC pc, PetscErrorCode (*)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx);
499: PETSC_EXTERN PetscErrorCode PCPatchSetComputeOperatorInteriorFacets(PC, PetscErrorCode (*)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
500: PETSC_EXTERN PetscErrorCode PCPatchSetComputeFunctionInteriorFacets(PC pc, PetscErrorCode (*)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx);

502: PETSC_EXTERN PetscErrorCode PCLMVMSetMatLMVM(PC, Mat);
503: PETSC_EXTERN PetscErrorCode PCLMVMGetMatLMVM(PC, Mat *);
504: PETSC_EXTERN PetscErrorCode PCLMVMSetIS(PC, IS);
505: PETSC_EXTERN PetscErrorCode PCLMVMClearIS(PC);
506: PETSC_EXTERN PetscErrorCode PCLMVMSetUpdateVec(PC, Vec);

508: PETSC_EXTERN PetscErrorCode PCExoticSetType(PC, PCExoticType);

510: PETSC_EXTERN PetscErrorCode PCDeflationSetInitOnly(PC, PetscBool);
511: PETSC_EXTERN PetscErrorCode PCDeflationSetLevels(PC, PetscInt);
512: PETSC_EXTERN PetscErrorCode PCDeflationSetReductionFactor(PC, PetscInt);
513: PETSC_EXTERN PetscErrorCode PCDeflationSetCorrectionFactor(PC, PetscScalar);
514: PETSC_EXTERN PetscErrorCode PCDeflationSetSpaceToCompute(PC, PCDeflationSpaceType, PetscInt);
515: PETSC_EXTERN PetscErrorCode PCDeflationSetSpace(PC, Mat, PetscBool);
516: PETSC_EXTERN PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC, Mat);
517: PETSC_EXTERN PetscErrorCode PCDeflationSetCoarseMat(PC, Mat);
518: PETSC_EXTERN PetscErrorCode PCDeflationGetPC(PC, PC *);

520: PETSC_EXTERN PetscErrorCode PCHPDDMSetAuxiliaryMat(PC, IS, Mat, PetscErrorCode (*)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void *), void *);
521: PETSC_EXTERN PetscErrorCode PCHPDDMSetRHSMat(PC, Mat);
522: PETSC_EXTERN PetscErrorCode PCHPDDMHasNeumannMat(PC, PetscBool);
523: PETSC_EXTERN PetscErrorCode PCHPDDMSetCoarseCorrectionType(PC, PCHPDDMCoarseCorrectionType);
524: PETSC_EXTERN PetscErrorCode PCHPDDMGetCoarseCorrectionType(PC, PCHPDDMCoarseCorrectionType *);
525: PETSC_EXTERN PetscErrorCode PCHPDDMSetSTShareSubKSP(PC, PetscBool);
526: PETSC_EXTERN PetscErrorCode PCHPDDMGetSTShareSubKSP(PC, PetscBool *);
527: PETSC_EXTERN PetscErrorCode PCHPDDMSetDeflationMat(PC, IS, Mat);
528: PETSC_EXTERN PetscErrorCode PCHPDDMFinalizePackage(void);
529: PETSC_EXTERN PetscErrorCode PCHPDDMInitializePackage(void);
530: PETSC_EXTERN PetscErrorCode PCHPDDMGetComplexities(PC, PetscReal *, PetscReal *);

532: PETSC_EXTERN PetscErrorCode PCAmgXGetResources(PC, void *);

534: PETSC_EXTERN PetscErrorCode PCMatSetApplyOperation(PC, MatOperation);
535: PETSC_EXTERN PetscErrorCode PCMatGetApplyOperation(PC, MatOperation *);