Actual source code: petscpc.h

petsc-3.6.1 2015-08-06
Report Typos and Errors
  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,PetscInt*);
 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 PCSetInitialGuessNonzero(PC,PetscBool);
 42: PETSC_EXTERN PetscErrorCode PCGetInitialGuessNonzero(PC,PetscBool*);
 43: PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool);
 44: PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*);


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

 49: PETSC_EXTERN PetscErrorCode PCReset(PC);
 50: PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
 51: PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);

 53: PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
 54: PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
 55: PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);

 57: PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat);
 58: PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*);
 59: PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);

 61: PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer);
 62: PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer);
 63: PETSC_STATIC_INLINE PetscErrorCode PCViewFromOptions(PC A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}

 65: PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
 66: PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
 67: PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);

 69: PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);

 71: /*
 72:       These are used to provide extra scaling of preconditioned
 73:    operator for time-stepping schemes like in SUNDIALS
 74: */
 75: PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *);
 76: PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
 77: PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
 78: PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec);

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

 82: PETSC_EXTERN PetscErrorCode PCJacobiSetType(PC,PCJacobiType);
 83: PETSC_EXTERN PetscErrorCode PCJacobiGetType(PC,PCJacobiType*);
 84: PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC,PetscBool);
 85: PETSC_EXTERN PetscErrorCode PCJacobiGetUseAbs(PC,PetscBool*);
 86: PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
 87: PETSC_EXTERN PetscErrorCode PCSORGetSymmetric(PC,MatSORType*);
 88: PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
 89: PETSC_EXTERN PetscErrorCode PCSORGetOmega(PC,PetscReal*);
 90: PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
 91: PETSC_EXTERN PetscErrorCode PCSORGetIterations(PC,PetscInt*,PetscInt*);

 93: PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
 94: PETSC_EXTERN PetscErrorCode PCEisenstatGetOmega(PC,PetscReal*);
 95: PETSC_EXTERN PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC,PetscBool);
 96: PETSC_EXTERN PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC,PetscBool*);

 98: PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
 99: PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
100: PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
101: PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);

103: PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
104: PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
105: PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
106: PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
107: PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
108: PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
109: PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
110: PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*);
111: PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**);
112: PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
113: PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]);

115: PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);

117: PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType);
118: PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal);

120: PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
121: PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
122: PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC);

124: PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
125: PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal);
126: PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);

128: PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType);
129: PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool );
130: PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool );
131: PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC,PetscBool);
132: PETSC_EXTERN PetscErrorCode PCFactorGetUseInPlace(PC,PetscBool*);
133: PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC,PetscBool);
134: PETSC_EXTERN PetscErrorCode PCFactorGetAllowDiagonalFill(PC,PetscBool*);
135: PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool);

137: PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
138: PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC,PetscInt*);
139: PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);

141: PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
142: PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
143: PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
144: PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool);
145: PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*);
146: PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool);

148: PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
149: PETSC_EXTERN PetscErrorCode PCASMGetType(PC,PCASMType*);
150: PETSC_EXTERN PetscErrorCode PCASMSetLocalType(PC,PCCompositeType);
151: PETSC_EXTERN PetscErrorCode PCASMGetLocalType(PC,PCCompositeType*);
152: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
153: PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
154: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
155: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
156: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);

158: PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt);
159: PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
160: PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
161: PETSC_EXTERN PetscErrorCode PCGASMSetUseDMSubdomains(PC,PetscBool);
162: PETSC_EXTERN PetscErrorCode PCGASMGetUseDMSubdomains(PC,PetscBool*);
163: PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );

165: PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
166: PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains(Mat,PetscInt,PetscInt*,IS*[]);
167: PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS*[],IS*[]);
168: PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
169: PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
170: PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);

172: PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
173: PETSC_EXTERN PetscErrorCode PCCompositeGetType(PC,PCCompositeType*);
174: PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
175: PETSC_EXTERN PetscErrorCode PCCompositeGetNumberPC(PC,PetscInt *);
176: PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
177: PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);

179: PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
180: PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
181: PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);

183: PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
184: PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
185: PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
186: PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
187: PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
188: PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
189: PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
190: PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);

192: PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
193: PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
194: PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteGradient(PC,Mat);
195: PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteCurl(PC,Mat);
196: PETSC_EXTERN PetscErrorCode PCHYPRESetEdgeConstantVectors(PC,Vec,Vec,Vec);
197: PETSC_EXTERN PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC,Mat);
198: PETSC_EXTERN PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC,Mat);
199: PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
200: PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);

202: PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
203: PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
204: PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
205: PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
206: PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
207: PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
208: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool);
209: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*);
210: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC,PetscBool);
211: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC,PetscBool*);
212: PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC,PetscBool);
213: PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC,PetscBool*);

215: PETSC_EXTERN PETSC_DEPRECATED("Use PCFieldSplitSetSchurPre") PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
216: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC,PCFieldSplitSchurPreType,Mat);
217: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC,PCFieldSplitSchurPreType*,Mat*);
218: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
219: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
220: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC,Mat *S);
221: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC,Mat *S);

223: PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
224: PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);

226: PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);

228: PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);

230: PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
231: PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);

233: PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
234: PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);

236: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
237: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
238: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);

240: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
241: PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
242: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
243: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);

245: PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC,PCPARMSGlobalType);
246: PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC,PCPARMSLocalType);
247: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC,PetscReal,PetscInt);
248: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC,PetscInt);
249: PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC,PetscBool);
250: PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC,PetscInt,PetscInt,PetscInt);

252: PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType);
253: PETSC_EXTERN PetscErrorCode PCGAMGGetType( PC,PCGAMGType*);
254: PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
255: PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
256: PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
257: PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
258: PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
259: PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
260: PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
261: PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC,PetscInt);
262: PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC,PetscBool);
263: PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscInt);
264: PETSC_EXTERN PetscErrorCode PCGAMGSetReuseInterpolation(PC,PetscBool);
265: PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void);
266: PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void);
267: PETSC_EXTERN PetscErrorCode PCGAMGRegister(PCGAMGType,PetscErrorCode (*)(PC));

269: PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC,PCGAMGClassicalType);
270: PETSC_EXTERN PetscErrorCode PCGAMGClassicalGetType(PC,PCGAMGClassicalType*);

272: PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisMat(PC,Mat);
273: PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC,IS);
274: PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
275: PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC,PetscInt);
276: PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace);
277: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
278: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC,IS);
279: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
280: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC,IS*);
281: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
282: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC,IS);
283: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
284: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC,IS*);
285: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
286: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC,PetscInt,IS[]);
287: PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode);
288: PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*);
289: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec);
290: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec);

292: PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool);
293: PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
294: PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec);

296: PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType);
297: PETSC_EXTERN PetscErrorCode PCMGGetType(PC,PCMGType*);
298: PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*);
299: PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*);

301: PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt);
302: PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt);
303: PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType);
304: PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType);
305: PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt);
306: PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt);
307: PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool);
308: PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*);

310: PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec);
311: PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec);
312: PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec);

314: PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat);
315: PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*);
316: PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat);
317: PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*);
318: PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec);
319: PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*);
320: PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
321: PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat,Vec,Vec,Vec);

323: #endif /* __PETSCPC_H */