Actual source code: petscpc.h
petsc-3.10.5 2019-03-28
1: /*
2: Preconditioner module.
3: */
6: #include <petscmat.h>
7: #include <petscpctypes.h>
9: PETSC_EXTERN PetscErrorCode PCInitializePackage(void);
11: /*
12: PCList contains the list of preconditioners currently registered
13: These are added with PCRegister()
14: */
15: PETSC_EXTERN PetscFunctionList PCList;
17: /* Logging support */
18: PETSC_EXTERN PetscClassId PC_CLASSID;
20: PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
21: PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType);
22: PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*);
23: PETSC_EXTERN PetscErrorCode PCSetUp(PC);
24: PETSC_EXTERN PetscErrorCode PCGetSetUpFailedReason(PC,PCFailedReason*);
25: PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
26: PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
27: PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
28: PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
29: PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
30: PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
31: PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *);
32: PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
33: PETSC_EXTERN PetscErrorCode PCSetReusePreconditioner(PC,PetscBool);
34: PETSC_EXTERN PetscErrorCode PCGetReusePreconditioner(PC,PetscBool*);
35: PETSC_EXTERN PetscErrorCode PCSetErrorIfFailure(PC,PetscBool);
37: #define PC_FILE_CLASSID 1211222
39: PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
40: PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *);
41: PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool);
42: PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*);
45: PETSC_EXTERN PetscErrorCode PCRegister(const char[],PetscErrorCode(*)(PC));
47: PETSC_EXTERN PetscErrorCode PCReset(PC);
48: PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
49: PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);
51: PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
52: PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
53: PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
55: PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat);
56: PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*);
57: PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
59: PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer);
60: PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer);
61: PETSC_STATIC_INLINE PetscErrorCode PCViewFromOptions(PC A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
63: PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
64: PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
65: PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);
67: PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
69: /*
70: These are used to provide extra scaling of preconditioned
71: operator for time-stepping schemes like in SUNDIALS
72: */
73: PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *);
74: PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
75: PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
76: PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec);
78: PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
79: PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);
81: PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
83: PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
84: PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);
86: /* ------------- options specific to particular preconditioners --------- */
88: PETSC_EXTERN PetscErrorCode PCJacobiSetType(PC,PCJacobiType);
89: PETSC_EXTERN PetscErrorCode PCJacobiGetType(PC,PCJacobiType*);
90: PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC,PetscBool);
91: PETSC_EXTERN PetscErrorCode PCJacobiGetUseAbs(PC,PetscBool*);
92: PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
93: PETSC_EXTERN PetscErrorCode PCSORGetSymmetric(PC,MatSORType*);
94: PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
95: PETSC_EXTERN PetscErrorCode PCSORGetOmega(PC,PetscReal*);
96: PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
97: PETSC_EXTERN PetscErrorCode PCSORGetIterations(PC,PetscInt*,PetscInt*);
99: PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
100: PETSC_EXTERN PetscErrorCode PCEisenstatGetOmega(PC,PetscReal*);
101: PETSC_EXTERN PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC,PetscBool);
102: PETSC_EXTERN PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC,PetscBool*);
104: PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
105: PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
106: PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
107: PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
109: PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
110: PETSC_EXTERN PetscErrorCode PCShellSetApplySymmetricLeft(PC,PetscErrorCode (*)(PC,Vec,Vec));
111: PETSC_EXTERN PetscErrorCode PCShellSetApplySymmetricRight(PC,PetscErrorCode (*)(PC,Vec,Vec));
112: PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
113: PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
114: PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
115: PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
116: PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
117: PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
118: PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*);
119: PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**);
120: PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
121: PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]);
123: PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);
125: PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType);
126: PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal);
128: PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverType(PC,MatSolverType);
129: PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverType(PC,MatSolverType*);
130: PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverType(PC);
131: PETSC_DEPRECATED("Use PCFactorSetMatSolverType") PETSC_STATIC_INLINE PetscErrorCode PCFactorSetMatSolverPackage(PC pc,MatSolverType stype)
132: { return PCFactorSetMatSolverType(pc,stype); }
133: PETSC_DEPRECATED("Use PCFactorGetMatSolverType") PETSC_STATIC_INLINE PetscErrorCode PCFactorGetMatSolverPackage(PC pc,MatSolverType *stype)
134: { return PCFactorGetMatSolverType(pc,stype); }
135: PETSC_DEPRECATED("Use PCFactorSetUpMatSolverType") PETSC_STATIC_INLINE PetscErrorCode PCFactorSetUpMatSolverPackage(PC pc)
136: { return PCFactorSetUpMatSolverType(pc); }
138: PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
139: PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal);
140: PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
142: PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType);
143: PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool );
144: PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool );
145: PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC,PetscBool);
146: PETSC_EXTERN PetscErrorCode PCFactorGetUseInPlace(PC,PetscBool*);
147: PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC,PetscBool);
148: PETSC_EXTERN PetscErrorCode PCFactorGetAllowDiagonalFill(PC,PetscBool*);
149: PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool);
151: PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
152: PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC,PetscInt*);
153: PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
154: PETSC_EXTERN PetscErrorCode PCFactorGetZeroPivot(PC,PetscReal*);
155: PETSC_EXTERN PetscErrorCode PCFactorGetShiftAmount(PC,PetscReal*);
156: PETSC_EXTERN PetscErrorCode PCFactorGetShiftType(PC,MatFactorShiftType*);
158: PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
159: PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
160: PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
161: PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool);
162: PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*);
163: PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool);
165: PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
166: PETSC_EXTERN PetscErrorCode PCASMGetType(PC,PCASMType*);
167: PETSC_EXTERN PetscErrorCode PCASMSetLocalType(PC,PCCompositeType);
168: PETSC_EXTERN PetscErrorCode PCASMGetLocalType(PC,PCCompositeType*);
169: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
170: PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
171: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
172: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
173: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
174: PETSC_EXTERN PetscErrorCode PCASMGetSubMatType(PC,MatType*);
175: PETSC_EXTERN PetscErrorCode PCASMSetSubMatType(PC,MatType);
177: PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt);
178: PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
179: PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
180: PETSC_EXTERN PetscErrorCode PCGASMSetUseDMSubdomains(PC,PetscBool);
181: PETSC_EXTERN PetscErrorCode PCGASMGetUseDMSubdomains(PC,PetscBool*);
182: PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );
184: PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
185: PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains(Mat,PetscInt,PetscInt*,IS*[]);
186: PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS*[],IS*[]);
187: PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
188: PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
189: PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);
191: PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
192: PETSC_EXTERN PetscErrorCode PCCompositeGetType(PC,PCCompositeType*);
193: PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
194: PETSC_EXTERN PetscErrorCode PCCompositeGetNumberPC(PC,PetscInt *);
195: PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
196: PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
198: PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
199: PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
200: PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
202: PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
203: PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
204: PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
205: PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
206: PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
207: PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
208: PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
209: PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
211: PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
212: PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
213: PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteGradient(PC,Mat);
214: PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteCurl(PC,Mat);
215: PETSC_EXTERN PetscErrorCode PCHYPRESetInterpolations(PC,PetscInt,Mat,Mat[],Mat,Mat[]);
216: PETSC_EXTERN PetscErrorCode PCHYPRESetEdgeConstantVectors(PC,Vec,Vec,Vec);
217: PETSC_EXTERN PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC,Mat);
218: PETSC_EXTERN PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC,Mat);
220: PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
221: PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
222: PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
223: PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
224: PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
225: PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
226: PETSC_EXTERN PetscErrorCode PCFieldSplitRestrictIS(PC,IS);
227: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool);
228: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*);
229: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC,PetscBool);
230: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC,PetscBool*);
231: PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC,PetscBool);
232: PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC,PetscBool*);
234: PETSC_EXTERN PETSC_DEPRECATED("Use PCFieldSplitSetSchurPre") PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
235: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC,PCFieldSplitSchurPreType,Mat);
236: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC,PCFieldSplitSchurPreType*,Mat*);
237: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
238: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurScale(PC,PetscScalar);
239: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
240: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC,Mat *S);
241: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC,Mat *S);
242: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDetectSaddlePoint(PC,PetscBool*);
243: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDetectSaddlePoint(PC,PetscBool);
245: PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
246: PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
247: PETSC_EXTERN PetscErrorCode PCGalerkinSetComputeSubmatrix(PC,PetscErrorCode (*)(PC,Mat,Mat,Mat*,void*),void*);
249: PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);
251: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
252: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
253: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
255: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
256: PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
257: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
258: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
260: PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC,PCPARMSGlobalType);
261: PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC,PCPARMSLocalType);
262: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC,PetscReal,PetscInt);
263: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC,PetscInt);
264: PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC,PetscBool);
265: PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC,PetscInt,PetscInt,PetscInt);
267: PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType);
268: PETSC_EXTERN PetscErrorCode PCGAMGGetType( PC,PCGAMGType*);
269: PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
270: PETSC_EXTERN PetscErrorCode PCGAMGSetRepartition(PC,PetscBool);
271: PETSC_EXTERN PetscErrorCode PCGAMGASMSetUseAggs(PC,PetscBool);
272: PETSC_EXTERN PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC,PetscBool);
273: PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
274: PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal[],PetscInt);
275: PETSC_EXTERN PetscErrorCode PCGAMGSetThresholdScale(PC,PetscReal);
276: PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
277: PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
278: PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC,PetscInt);
279: PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC,PetscBool);
280: PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscInt);
281: PETSC_EXTERN PetscErrorCode PCGAMGSetReuseInterpolation(PC,PetscBool);
282: PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void);
283: PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void);
284: PETSC_EXTERN PetscErrorCode PCGAMGRegister(PCGAMGType,PetscErrorCode (*)(PC));
286: PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC,PCGAMGClassicalType);
287: PETSC_EXTERN PetscErrorCode PCGAMGClassicalGetType(PC,PCGAMGClassicalType*);
289: PETSC_EXTERN PetscErrorCode PCBDDCSetDiscreteGradient(PC,Mat,PetscInt,PetscInt,PetscBool,PetscBool);
290: PETSC_EXTERN PetscErrorCode PCBDDCSetDivergenceMat(PC,Mat,PetscBool,IS);
291: PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisMat(PC,Mat,PetscBool);
292: PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesIS(PC,IS);
293: PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC,IS);
294: PETSC_EXTERN PetscErrorCode PCBDDCGetPrimalVerticesIS(PC,IS*);
295: PETSC_EXTERN PetscErrorCode PCBDDCGetPrimalVerticesLocalIS(PC,IS*);
296: PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
297: PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC,PetscInt);
298: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
299: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC,IS);
300: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
301: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC,IS*);
302: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
303: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC,IS);
304: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
305: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC,IS*);
306: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
307: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC,PetscInt,IS[]);
308: PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode);
309: PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,PetscBool,const char*,Mat*,PC*);
310: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec);
311: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec);
312: PETSC_EXTERN PetscErrorCode PCBDDCFinalizePackage(void);
313: PETSC_EXTERN PetscErrorCode PCBDDCInitializePackage(void);
315: PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool);
316: PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
317: PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec);
319: PETSC_EXTERN PetscInt PetscMGLevelId;
320: PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType);
321: PETSC_EXTERN PetscErrorCode PCMGGetType(PC,PCMGType*);
322: PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*);
323: PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*);
325: PETSC_EXTERN PetscErrorCode PCMGSetDistinctSmoothUp(PC);
326: PETSC_EXTERN PetscErrorCode PCMGSetNumberSmooth(PC,PetscInt);
327: PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType);
328: PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType);
329: PETSC_DEPRECATED("Use PCMGSetCycleTypeOnLevel") PETSC_STATIC_INLINE PetscErrorCode PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt t) {return PCMGSetCycleTypeOnLevel(pc,l,(PCMGCycleType)t);}
330: PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt);
331: PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PCMGGalerkinType);
332: PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PCMGGalerkinType*);
334: PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec);
335: PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec);
336: PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec);
338: PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat);
339: PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*);
340: PETSC_EXTERN PetscErrorCode PCMGSetInjection(PC,PetscInt,Mat);
341: PETSC_EXTERN PetscErrorCode PCMGGetInjection(PC,PetscInt,Mat*);
342: PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat);
343: PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*);
344: PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec);
345: PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*);
346: PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
347: PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat,Vec,Vec,Vec);
349: PETSC_EXTERN PetscErrorCode PCTelescopeGetSubcommType(PC,PetscSubcommType*);
350: PETSC_EXTERN PetscErrorCode PCTelescopeSetSubcommType(PC,PetscSubcommType);
351: PETSC_EXTERN PetscErrorCode PCTelescopeGetReductionFactor(PC,PetscInt*);
352: PETSC_EXTERN PetscErrorCode PCTelescopeSetReductionFactor(PC,PetscInt);
353: PETSC_EXTERN PetscErrorCode PCTelescopeGetIgnoreDM(PC,PetscBool*);
354: PETSC_EXTERN PetscErrorCode PCTelescopeSetIgnoreDM(PC,PetscBool);
355: PETSC_EXTERN PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC,PetscBool*);
356: PETSC_EXTERN PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC,PetscBool);
357: PETSC_EXTERN PetscErrorCode PCTelescopeGetDM(PC,DM*);
359: PETSC_EXTERN PetscErrorCode PCPatchSetSaveOperators(PC, PetscBool);
360: PETSC_EXTERN PetscErrorCode PCPatchGetSaveOperators(PC, PetscBool *);
361: PETSC_EXTERN PetscErrorCode PCPatchSetPartitionOfUnity(PC, PetscBool);
362: PETSC_EXTERN PetscErrorCode PCPatchGetPartitionOfUnity(PC, PetscBool *);
363: PETSC_EXTERN PetscErrorCode PCPatchSetMultiplicative(PC, PetscBool);
364: PETSC_EXTERN PetscErrorCode PCPatchGetMultiplicative(PC, PetscBool *);
365: PETSC_EXTERN PetscErrorCode PCPatchSetSubMatType(PC, MatType);
366: PETSC_EXTERN PetscErrorCode PCPatchGetSubMatType(PC, MatType *);
367: PETSC_EXTERN PetscErrorCode PCPatchSetCellNumbering(PC, PetscSection);
368: PETSC_EXTERN PetscErrorCode PCPatchGetCellNumbering(PC, PetscSection *);
369: PETSC_EXTERN PetscErrorCode PCPatchSetConstructType(PC, PCPatchConstructType, PetscErrorCode (*)(PC, PetscInt *, IS **, IS *, void *), void *);
370: PETSC_EXTERN PetscErrorCode PCPatchGetConstructType(PC, PCPatchConstructType *, PetscErrorCode (**)(PC, PetscInt *, IS **, IS *, void *), void **);
371: PETSC_EXTERN PetscErrorCode PCPatchSetDiscretisationInfo(PC, PetscInt, DM *, PetscInt *, PetscInt *, const PetscInt **, const PetscInt *, PetscInt, const PetscInt *, PetscInt, const PetscInt *);
372: PETSC_EXTERN PetscErrorCode PCPatchSetComputeOperator(PC, PetscErrorCode (*)(PC,PetscInt,Mat,IS,PetscInt,const PetscInt *,void *), void *);
374: PETSC_EXTERN PetscErrorCode PCLMVMSetMatLMVM(PC, Mat);
375: PETSC_EXTERN PetscErrorCode PCLMVMGetMatLMVM(PC, Mat*);
376: PETSC_EXTERN PetscErrorCode PCLMVMSetIS(PC, IS);
377: PETSC_EXTERN PetscErrorCode PCLMVMClearIS(PC);
379: #endif /* __PETSCPC_H */