Actual source code: petscpc.h

petsc-3.12.5 2020-03-29
Report Typos and Errors
  1: /*
  2:       Preconditioner module.
  3: */
  4: #if !defined(PETSCPC_H)
  5: #define PETSCPC_H
  6:  #include <petscmat.h>
  7:  #include <petscdmtypes.h>
  8:  #include <petscpctypes.h>

 10: PETSC_EXTERN PetscErrorCode PCInitializePackage(void);

 12: /*
 13:     PCList contains the list of preconditioners currently registered
 14:    These are added with PCRegister()
 15: */
 16: PETSC_EXTERN PetscFunctionList PCList;

 18: /* Logging support */
 19: PETSC_EXTERN PetscClassId PC_CLASSID;

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

 39: PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
 40: PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType);
 41: PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*);
 42: PETSC_EXTERN PetscErrorCode PCSetUp(PC);
 43: PETSC_EXTERN PetscErrorCode PCGetFailedReason(PC,PCFailedReason*);
 44: PETSC_DEPRECATED_FUNCTION("Use PCGetFailedReason() (since version 3.11)") PETSC_STATIC_INLINE PetscErrorCode PCGetSetUpFailedReason(PC pc,PCFailedReason *reason) {return PCGetFailedReason(pc,reason);}
 45: PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
 46: PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
 47: PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
 48: PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
 49: PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
 50: PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
 51: PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *);
 52: PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
 53: PETSC_EXTERN PetscErrorCode PCSetReusePreconditioner(PC,PetscBool);
 54: PETSC_EXTERN PetscErrorCode PCGetReusePreconditioner(PC,PetscBool*);
 55: PETSC_EXTERN PetscErrorCode PCSetErrorIfFailure(PC,PetscBool);

 57: #define PC_FILE_CLASSID 1211222

 59: PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
 60: PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *);
 61: PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool);
 62: PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*);


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

 67: PETSC_EXTERN PetscErrorCode PCReset(PC);
 68: PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
 69: PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);

 71: PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
 72: PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
 73: PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);

 75: PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat);
 76: PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*);
 77: PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);

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

 83: PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
 84: PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
 85: PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);

 87: PETSC_EXTERN PetscErrorCode PCComputeOperator(PC,MatType,Mat*);
 88: PETSC_DEPRECATED_FUNCTION("Use PCComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode PCComputeExplicitOperator(PC A,Mat* B) { return PCComputeOperator(A,NULL,B); }

 90: /*
 91:       These are used to provide extra scaling of preconditioned
 92:    operator for time-stepping schemes like in SUNDIALS
 93: */
 94: PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *);
 95: PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
 96: PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
 97: PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec);

 99: PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
100: PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);

102: PETSC_EXTERN PetscErrorCode PCGetInterpolations(PC,PetscInt*,Mat*[]);
103: PETSC_EXTERN PetscErrorCode PCGetCoarseOperators(PC pc,PetscInt*,Mat*[]);

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

107: PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
108: PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);

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

112: PETSC_EXTERN PetscErrorCode PCJacobiSetType(PC,PCJacobiType);
113: PETSC_EXTERN PetscErrorCode PCJacobiGetType(PC,PCJacobiType*);
114: PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC,PetscBool);
115: PETSC_EXTERN PetscErrorCode PCJacobiGetUseAbs(PC,PetscBool*);
116: PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
117: PETSC_EXTERN PetscErrorCode PCSORGetSymmetric(PC,MatSORType*);
118: PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
119: PETSC_EXTERN PetscErrorCode PCSORGetOmega(PC,PetscReal*);
120: PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
121: PETSC_EXTERN PetscErrorCode PCSORGetIterations(PC,PetscInt*,PetscInt*);

123: PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
124: PETSC_EXTERN PetscErrorCode PCEisenstatGetOmega(PC,PetscReal*);
125: PETSC_EXTERN PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC,PetscBool);
126: PETSC_EXTERN PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC,PetscBool*);

128: PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
129: PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
130: PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
131: PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);

133: PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
134: PETSC_EXTERN PetscErrorCode PCShellSetApplySymmetricLeft(PC,PetscErrorCode (*)(PC,Vec,Vec));
135: PETSC_EXTERN PetscErrorCode PCShellSetApplySymmetricRight(PC,PetscErrorCode (*)(PC,Vec,Vec));
136: PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
137: PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
138: PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
139: PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
140: PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
141: PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
142: PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*);
143: PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**);
144: PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
145: PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]);

147: PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);

149: PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType);
150: PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal);

152: PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverType(PC,MatSolverType);
153: PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverType(PC,MatSolverType*);
154: PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverType(PC);
155: PETSC_DEPRECATED_FUNCTION("Use PCFactorSetMatSolverType() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode PCFactorSetMatSolverPackage(PC pc,MatSolverType stype) {return PCFactorSetMatSolverType(pc,stype);}
156: PETSC_DEPRECATED_FUNCTION("Use PCFactorGetMatSolverType() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode PCFactorGetMatSolverPackage(PC pc,MatSolverType *stype) {return PCFactorGetMatSolverType(pc,stype);}
157: PETSC_DEPRECATED_FUNCTION("Use PCFactorSetUpMatSolverType() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode PCFactorSetUpMatSolverPackage(PC pc) {return PCFactorSetUpMatSolverType(pc);}

159: PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
160: PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal);
161: PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);

163: PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType);
164: PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool );
165: PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool );
166: PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC,PetscBool);
167: PETSC_EXTERN PetscErrorCode PCFactorGetUseInPlace(PC,PetscBool*);
168: PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC,PetscBool);
169: PETSC_EXTERN PetscErrorCode PCFactorGetAllowDiagonalFill(PC,PetscBool*);
170: PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool);

172: PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
173: PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC,PetscInt*);
174: PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
175: PETSC_EXTERN PetscErrorCode PCFactorGetZeroPivot(PC,PetscReal*);
176: PETSC_EXTERN PetscErrorCode PCFactorGetShiftAmount(PC,PetscReal*);
177: PETSC_EXTERN PetscErrorCode PCFactorGetShiftType(PC,MatFactorShiftType*);

179: PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
180: PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
181: PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
182: PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool);
183: PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*);
184: PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool);

186: PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
187: PETSC_EXTERN PetscErrorCode PCASMGetType(PC,PCASMType*);
188: PETSC_EXTERN PetscErrorCode PCASMSetLocalType(PC,PCCompositeType);
189: PETSC_EXTERN PetscErrorCode PCASMGetLocalType(PC,PCCompositeType*);
190: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
191: PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
192: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
193: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
194: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
195: PETSC_EXTERN PetscErrorCode PCASMGetSubMatType(PC,MatType*);
196: PETSC_EXTERN PetscErrorCode PCASMSetSubMatType(PC,MatType);

198: PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt);
199: PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
200: PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
201: PETSC_EXTERN PetscErrorCode PCGASMSetUseDMSubdomains(PC,PetscBool);
202: PETSC_EXTERN PetscErrorCode PCGASMGetUseDMSubdomains(PC,PetscBool*);
203: PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );

205: PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
206: PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains(Mat,PetscInt,PetscInt*,IS*[]);
207: PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS*[],IS*[]);
208: PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
209: PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
210: PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);

212: PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
213: PETSC_EXTERN PetscErrorCode PCCompositeGetType(PC,PCCompositeType*);
214: PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
215: PETSC_EXTERN PetscErrorCode PCCompositeGetNumberPC(PC,PetscInt *);
216: PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
217: PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);

219: PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
220: PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
221: PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);

223: PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
224: PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
225: PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
226: PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
227: PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
228: PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
229: PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
230: PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);

232: PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
233: PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
234: PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteGradient(PC,Mat);
235: PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteCurl(PC,Mat);
236: PETSC_EXTERN PetscErrorCode PCHYPRESetInterpolations(PC,PetscInt,Mat,Mat[],Mat,Mat[]);
237: PETSC_EXTERN PetscErrorCode PCHYPRESetEdgeConstantVectors(PC,Vec,Vec,Vec);
238: PETSC_EXTERN PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC,Mat);
239: PETSC_EXTERN PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC,Mat);

241: PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
242: PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
243: PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
244: PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
245: PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
246: PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
247: PETSC_EXTERN PetscErrorCode PCFieldSplitGetISByIndex(PC,PetscInt,IS*);
248: PETSC_EXTERN PetscErrorCode PCFieldSplitRestrictIS(PC,IS);
249: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool);
250: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*);
251: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC,PetscBool);
252: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC,PetscBool*);
253: PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC,PetscBool);
254: PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC,PetscBool*);

256: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use PCFieldSplitSetSchurPre() (since version 3.5)") PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
257: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC,PCFieldSplitSchurPreType,Mat);
258: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC,PCFieldSplitSchurPreType*,Mat*);
259: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
260: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurScale(PC,PetscScalar);
261: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
262: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC,Mat *S);
263: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC,Mat *S);
264: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDetectSaddlePoint(PC,PetscBool*);
265: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDetectSaddlePoint(PC,PetscBool);
266: PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBTol(PC,PetscReal);
267: PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBNu(PC,PetscReal);
268: PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBMaxit(PC,PetscInt);
269: PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBDelay(PC,PetscInt);

271: PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
272: PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
273: PETSC_EXTERN PetscErrorCode PCGalerkinSetComputeSubmatrix(PC,PetscErrorCode (*)(PC,Mat,Mat,Mat*,void*),void*);

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

277: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
278: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
279: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);

281: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
282: PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
283: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
284: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);

286: PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC,PCPARMSGlobalType);
287: PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC,PCPARMSLocalType);
288: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC,PetscReal,PetscInt);
289: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC,PetscInt);
290: PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC,PetscBool);
291: PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC,PetscInt,PetscInt,PetscInt);

293: PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType);
294: PETSC_EXTERN PetscErrorCode PCGAMGGetType( PC,PCGAMGType*);
295: PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
296: PETSC_EXTERN PetscErrorCode PCGAMGSetRepartition(PC,PetscBool);
297: PETSC_EXTERN PetscErrorCode PCGAMGASMSetUseAggs(PC,PetscBool);
298: PETSC_EXTERN PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC,PetscBool);
299: PETSC_EXTERN PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC,PetscBool);
300: PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC,PCGAMGLayoutType);
301: PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
302: PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal[],PetscInt);
303: PETSC_EXTERN PetscErrorCode PCGAMGSetThresholdScale(PC,PetscReal);
304: PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
305: PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
306: PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC,PetscInt);
307: PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC,PetscBool);
308: PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscInt);
309: PETSC_EXTERN PetscErrorCode PCGAMGSetReuseInterpolation(PC,PetscBool);
310: PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void);
311: PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void);
312: PETSC_EXTERN PetscErrorCode PCGAMGRegister(PCGAMGType,PetscErrorCode (*)(PC));

314: PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC,PCGAMGClassicalType);
315: PETSC_EXTERN PetscErrorCode PCGAMGClassicalGetType(PC,PCGAMGClassicalType*);

317: PETSC_EXTERN PetscErrorCode PCBDDCSetDiscreteGradient(PC,Mat,PetscInt,PetscInt,PetscBool,PetscBool);
318: PETSC_EXTERN PetscErrorCode PCBDDCSetDivergenceMat(PC,Mat,PetscBool,IS);
319: PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisMat(PC,Mat,PetscBool);
320: PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesIS(PC,IS);
321: PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC,IS);
322: PETSC_EXTERN PetscErrorCode PCBDDCGetPrimalVerticesIS(PC,IS*);
323: PETSC_EXTERN PetscErrorCode PCBDDCGetPrimalVerticesLocalIS(PC,IS*);
324: PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
325: PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC,PetscInt);
326: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
327: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC,IS);
328: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
329: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC,IS*);
330: PETSC_EXTERN PetscErrorCode PCBDDCSetInterfaceExtType(PC,PCBDDCInterfaceExtType);
331: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
332: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC,IS);
333: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
334: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC,IS*);
335: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
336: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC,PetscInt,IS[]);
337: PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode);
338: PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,PetscBool,const char*,Mat*,PC*);
339: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec);
340: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec);
341: PETSC_EXTERN PetscErrorCode PCBDDCFinalizePackage(void);
342: PETSC_EXTERN PetscErrorCode PCBDDCInitializePackage(void);

344: PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool);
345: PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
346: PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec);

348: PETSC_EXTERN PetscInt PetscMGLevelId;
349: PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType);
350: PETSC_EXTERN PetscErrorCode PCMGGetType(PC,PCMGType*);
351: PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*);
352: PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*);

354: PETSC_EXTERN PetscErrorCode PCMGSetDistinctSmoothUp(PC);
355: PETSC_EXTERN PetscErrorCode PCMGSetNumberSmooth(PC,PetscInt);
356: PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType);
357: PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType);
358: PETSC_DEPRECATED_FUNCTION("Use PCMGSetCycleTypeOnLevel() (since version 3.5)") PETSC_STATIC_INLINE PetscErrorCode PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt t) {return PCMGSetCycleTypeOnLevel(pc,l,(PCMGCycleType)t);}
359: PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt);
360: PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PCMGGalerkinType);
361: PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PCMGGalerkinType*);

363: PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec);
364: PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec);
365: PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec);

367: PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat);
368: PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*);
369: PETSC_EXTERN PetscErrorCode PCMGSetInjection(PC,PetscInt,Mat);
370: PETSC_EXTERN PetscErrorCode PCMGGetInjection(PC,PetscInt,Mat*);
371: PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat);
372: PETSC_EXTERN PetscErrorCode PCMGSetOperators(PC,PetscInt,Mat,Mat);
373: PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*);
374: PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec);
375: PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*);
376: PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
377: PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat,Vec,Vec,Vec);

379: PETSC_EXTERN PetscErrorCode PCHMGSetReuseInterpolation(PC,PetscBool);
380: PETSC_EXTERN PetscErrorCode PCHMGSetUseSubspaceCoarsening(PC,PetscBool);
381: PETSC_EXTERN PetscErrorCode PCHMGSetInnerPCType(PC,PCType);
382: PETSC_EXTERN PetscErrorCode PCHMGSetCoarseningComponent(PC,PetscInt);
383: PETSC_EXTERN PetscErrorCode PCHMGUseMatMAIJ(PC,PetscBool);

385: PETSC_EXTERN PetscErrorCode PCTelescopeGetSubcommType(PC,PetscSubcommType*);
386: PETSC_EXTERN PetscErrorCode PCTelescopeSetSubcommType(PC,PetscSubcommType);
387: PETSC_EXTERN PetscErrorCode PCTelescopeGetReductionFactor(PC,PetscInt*);
388: PETSC_EXTERN PetscErrorCode PCTelescopeSetReductionFactor(PC,PetscInt);
389: PETSC_EXTERN PetscErrorCode PCTelescopeGetIgnoreDM(PC,PetscBool*);
390: PETSC_EXTERN PetscErrorCode PCTelescopeSetIgnoreDM(PC,PetscBool);
391: PETSC_EXTERN PetscErrorCode PCTelescopeGetUseCoarseDM(PC,PetscBool*);
392: PETSC_EXTERN PetscErrorCode PCTelescopeSetUseCoarseDM(PC,PetscBool);
393: PETSC_EXTERN PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC,PetscBool*);
394: PETSC_EXTERN PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC,PetscBool);
395: PETSC_EXTERN PetscErrorCode PCTelescopeGetDM(PC,DM*);

397: PETSC_EXTERN PetscErrorCode PCPatchSetSaveOperators(PC, PetscBool);
398: PETSC_EXTERN PetscErrorCode PCPatchGetSaveOperators(PC, PetscBool *);
399: PETSC_EXTERN PetscErrorCode PCPatchSetPrecomputeElementTensors(PC, PetscBool);
400: PETSC_EXTERN PetscErrorCode PCPatchGetPrecomputeElementTensors(PC, PetscBool *);
401: PETSC_EXTERN PetscErrorCode PCPatchSetPartitionOfUnity(PC, PetscBool);
402: PETSC_EXTERN PetscErrorCode PCPatchGetPartitionOfUnity(PC, PetscBool *);
403: PETSC_EXTERN PetscErrorCode PCPatchSetMultiplicative(PC, PetscBool);
404: PETSC_EXTERN PetscErrorCode PCPatchGetMultiplicative(PC, PetscBool *);
405: PETSC_EXTERN PetscErrorCode PCPatchSetSubMatType(PC, MatType);
406: PETSC_EXTERN PetscErrorCode PCPatchGetSubMatType(PC, MatType *);
407: PETSC_EXTERN PetscErrorCode PCPatchSetCellNumbering(PC, PetscSection);
408: PETSC_EXTERN PetscErrorCode PCPatchGetCellNumbering(PC, PetscSection *);
409: PETSC_EXTERN PetscErrorCode PCPatchSetConstructType(PC, PCPatchConstructType,   PetscErrorCode (*)(PC, PetscInt *, IS **, IS *, void *), void *);
410: PETSC_EXTERN PetscErrorCode PCPatchGetConstructType(PC, PCPatchConstructType *, PetscErrorCode (**)(PC, PetscInt *, IS **, IS *, void *), void **);
411: PETSC_EXTERN PetscErrorCode PCPatchSetDiscretisationInfo(PC, PetscInt, DM *, PetscInt *, PetscInt *, const PetscInt **, const PetscInt *, PetscInt, const PetscInt *, PetscInt, const PetscInt *);
412: PETSC_EXTERN PetscErrorCode PCPatchSetComputeOperator(PC, PetscErrorCode (*)(PC,PetscInt,Vec,Mat,IS,PetscInt,const PetscInt *,const PetscInt *, void *), void *);
413: PETSC_EXTERN PetscErrorCode PCPatchSetComputeFunction(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx);
414: PETSC_EXTERN PetscErrorCode PCPatchSetComputeOperatorInteriorFacets(PC, PetscErrorCode (*)(PC,PetscInt,Vec,Mat,IS,PetscInt,const PetscInt *,const PetscInt *, void *), void *);
415: PETSC_EXTERN PetscErrorCode PCPatchSetComputeFunctionInteriorFacets(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx);

417: PETSC_EXTERN PetscErrorCode PCLMVMSetMatLMVM(PC, Mat);
418: PETSC_EXTERN PetscErrorCode PCLMVMGetMatLMVM(PC, Mat*);
419: PETSC_EXTERN PetscErrorCode PCLMVMSetIS(PC, IS);
420: PETSC_EXTERN PetscErrorCode PCLMVMClearIS(PC);

422: PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);

424: PETSC_EXTERN PetscErrorCode PCDeflationSetInitOnly(PC,PetscBool);
425: PETSC_EXTERN PetscErrorCode PCDeflationSetLevels(PC,PetscInt);
426: PETSC_EXTERN PetscErrorCode PCDeflationSetReductionFactor(PC,PetscInt);
427: PETSC_EXTERN PetscErrorCode PCDeflationSetCorrectionFactor(PC,PetscScalar);
428: PETSC_EXTERN PetscErrorCode PCDeflationSetSpaceToCompute(PC,PCDeflationSpaceType,PetscInt);
429: PETSC_EXTERN PetscErrorCode PCDeflationSetSpace(PC,Mat,PetscBool);
430: PETSC_EXTERN PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC,Mat);
431: PETSC_EXTERN PetscErrorCode PCDeflationSetCoarseMat(PC,Mat);
432: PETSC_EXTERN PetscErrorCode PCDeflationGetPC(PC,PC*);

434: PETSC_EXTERN PetscErrorCode PCHPDDMSetAuxiliaryMat(PC,IS,Mat,PetscErrorCode (*)(Mat,PetscReal,Vec,Vec,PetscReal,IS,void*),void*);
435: PETSC_EXTERN PetscErrorCode PCHPDDMSetCoarseCorrectionType(PC,PCHPDDMCoarseCorrectionType);
436: PETSC_EXTERN PetscErrorCode PCHPDDMGetCoarseCorrectionType(PC,PCHPDDMCoarseCorrectionType*);

438: #endif /* PETSCPC_H */