Actual source code: petscpc.h
petsc-3.5.0 2014-06-30
1: /*
2: Preconditioner module.
3: */
6: #include <petscmat.h>
7: #include <petscdmtypes.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: /*S
18: PC - Abstract PETSc object that manages all preconditioners including direct solvers such as PCLU
20: Level: beginner
22: Concepts: preconditioners
24: .seealso: PCCreate(), PCSetType(), PCType (for list of available types)
25: S*/
26: typedef struct _p_PC* PC;
28: /*J
29: PCType - String with the name of a PETSc preconditioner method.
31: Level: beginner
33: Notes: Click on the links below to see details on a particular solver
35: PCRegister() is used to register preconditioners that are then accessible via PCSetType()
37: .seealso: PCSetType(), PC, PCCreate(), PCRegister(), PCSetFromOptions()
38: J*/
39: typedef const char* PCType;
40: #define PCNONE "none"
41: #define PCJACOBI "jacobi"
42: #define PCSOR "sor"
43: #define PCLU "lu"
44: #define PCSHELL "shell"
45: #define PCBJACOBI "bjacobi"
46: #define PCMG "mg"
47: #define PCEISENSTAT "eisenstat"
48: #define PCILU "ilu"
49: #define PCICC "icc"
50: #define PCASM "asm"
51: #define PCGASM "gasm"
52: #define PCKSP "ksp"
53: #define PCCOMPOSITE "composite"
54: #define PCREDUNDANT "redundant"
55: #define PCSPAI "spai"
56: #define PCNN "nn"
57: #define PCCHOLESKY "cholesky"
58: #define PCPBJACOBI "pbjacobi"
59: #define PCMAT "mat"
60: #define PCHYPRE "hypre"
61: #define PCPARMS "parms"
62: #define PCFIELDSPLIT "fieldsplit"
63: #define PCTFS "tfs"
64: #define PCML "ml"
65: #define PCGALERKIN "galerkin"
66: #define PCEXOTIC "exotic"
67: #define PCCP "cp"
68: #define PCBFBT "bfbt"
69: #define PCLSC "lsc"
70: #define PCPYTHON "python"
71: #define PCPFMG "pfmg"
72: #define PCSYSPFMG "syspfmg"
73: #define PCREDISTRIBUTE "redistribute"
74: #define PCSVD "svd"
75: #define PCGAMG "gamg"
76: #define PCSACUSP "sacusp" /* these four run on NVIDIA GPUs using CUSP */
77: #define PCSACUSPPOLY "sacusppoly"
78: #define PCBICGSTABCUSP "bicgstabcusp"
79: #define PCAINVCUSP "ainvcusp"
80: #define PCBDDC "bddc"
81: #define PCKACZMARZ "kaczmarz"
83: /* Logging support */
84: PETSC_EXTERN PetscClassId PC_CLASSID;
86: /*E
87: PCSide - If the preconditioner is to be applied to the left, right
88: or symmetrically around the operator.
90: Level: beginner
92: .seealso:
93: E*/
94: typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
95: #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
96: PETSC_EXTERN const char *const *const PCSides;
98: PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
99: PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType);
100: PETSC_EXTERN PetscErrorCode PCSetUp(PC);
101: PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
102: PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
103: PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
104: PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
105: PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
106: PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
107: PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *);
108: PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
109: PETSC_EXTERN PetscErrorCode PCSetReusePreconditioner(PC,PetscBool);
111: #define PC_FILE_CLASSID 1211222
113: /*E
114: PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
116: Level: advanced
118: Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
120: .seealso: PCApplyRichardson()
121: E*/
122: typedef enum {
123: PCRICHARDSON_CONVERGED_RTOL = 2,
124: PCRICHARDSON_CONVERGED_ATOL = 3,
125: PCRICHARDSON_CONVERGED_ITS = 4,
126: PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason;
128: PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
129: PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *);
130: PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool );
131: PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool);
132: PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*);
134: PETSC_EXTERN PetscErrorCode PCRegisterAll(void);
135: PETSC_EXTERN PetscBool PCRegisterAllCalled;
137: PETSC_EXTERN PetscErrorCode PCRegister(const char[],PetscErrorCode(*)(PC));
139: PETSC_EXTERN PetscErrorCode PCReset(PC);
140: PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
141: PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);
142: PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*);
144: PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
145: PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
146: PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
148: PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat);
149: PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*);
150: PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
152: PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer);
153: PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer);
154: PETSC_STATIC_INLINE PetscErrorCode PCViewFromOptions(PC A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
156: PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
157: PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
158: PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);
160: PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
162: /*
163: These are used to provide extra scaling of preconditioned
164: operator for time-stepping schemes like in SUNDIALS
165: */
166: PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *);
167: PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
168: PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
169: PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec);
171: /* ------------- options specific to particular preconditioners --------- */
173: PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
174: PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC);
175: PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC);
176: PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
177: PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
178: PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
180: PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
181: PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
183: PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
184: PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
186: PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
187: PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
188: PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
189: PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
190: PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
191: PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
192: PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
193: PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**);
194: PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*);
195: PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
196: PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]);
198: PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);
200: PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType);
201: PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal);
203: PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
204: PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
205: PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC);
207: PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
208: PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal);
209: PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
211: PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType);
212: PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool );
213: PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool );
214: PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC);
215: PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
216: PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool );
218: PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC,PetscInt*);
219: PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
220: PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
222: PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
223: PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
224: PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
225: PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool);
226: PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*);
227: PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool );
229: /*E
230: PCASMType - Type of additive Schwarz method to use
232: $ PC_ASM_BASIC - Symmetric version where residuals from the ghost points are used
233: $ and computed values in ghost regions are added together.
234: $ Classical standard additive Schwarz.
235: $ PC_ASM_RESTRICT - Residuals from ghost points are used but computed values in ghost
236: $ region are discarded.
237: $ Default.
238: $ PC_ASM_INTERPOLATE - Residuals from ghost points are not used, computed values in ghost
239: $ region are added back in.
240: $ PC_ASM_NONE - Residuals from ghost points are not used, computed ghost values are
241: $ discarded.
242: $ Not very good.
244: Level: beginner
246: .seealso: PCASMSetType()
247: E*/
248: typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
249: PETSC_EXTERN const char *const PCASMTypes[];
251: PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
252: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
253: PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
254: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
255: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
256: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
258: /*E
259: PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
261: Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational
262: domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute
263: a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
264: (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
266: $ PC_GASM_BASIC - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
267: $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections
268: $ from neighboring subdomains.
269: $ Classical standard additive Schwarz.
270: $ PC_GASM_RESTRICT - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
271: $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result,
272: $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
273: $ assumption).
274: $ Default.
275: $ PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
276: $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
277: $ from neighboring subdomains.
278: $
279: $ PC_GASM_NONE - Residuals and corrections are zeroed out outside the local subdomains.
280: $ Not very good.
282: Level: beginner
284: .seealso: PCGASMSetType()
285: E*/
286: typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
287: PETSC_EXTERN const char *const PCGASMTypes[];
289: PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
290: PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool);
291: PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
292: PETSC_EXTERN PetscErrorCode PCGASMSetDMSubdomains(PC,PetscBool);
293: PETSC_EXTERN PetscErrorCode PCGASMGetDMSubdomains(PC,PetscBool*);
294: PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );
296: PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
297: PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]);
298: PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
299: PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
300: PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
301: PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);
303: /*E
304: PCCompositeType - Determines how two or more preconditioner are composed
306: $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
307: $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
308: $ computed after the previous preconditioner application
309: $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
310: $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
311: $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
312: $ where first preconditioner is built from alpha I + S and second from
313: $ alpha I + R
315: Level: beginner
317: .seealso: PCCompositeSetType()
318: E*/
319: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
320: PETSC_EXTERN const char *const PCCompositeTypes[];
322: PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
323: PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
324: PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
325: PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
327: PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
328: PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
329: PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
331: PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
332: PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
333: PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
334: PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
335: PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
336: PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
337: PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
338: PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
340: PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
341: PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
342: PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
343: PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
345: PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
346: PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
347: PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
348: PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
349: PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
350: PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
351: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool);
352: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*);
353: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC,PetscBool);
354: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC,PetscBool*);
355: PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC,PetscBool);
356: PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC,PetscBool*);
359: /*E
360: PCFieldSplitSchurPreType - Determines how to precondition Schur complement
362: Level: intermediate
364: .seealso: PCFieldSplitSetSchurPre()
365: E*/
366: typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_SELFP,PC_FIELDSPLIT_SCHUR_PRE_A11,PC_FIELDSPLIT_SCHUR_PRE_USER,PC_FIELDSPLIT_SCHUR_PRE_FULL} PCFieldSplitSchurPreType;
367: PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
369: /*E
370: PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
372: Level: intermediate
374: .seealso: PCFieldSplitSetSchurFactType()
375: E*/
376: typedef enum {
377: PC_FIELDSPLIT_SCHUR_FACT_DIAG,
378: PC_FIELDSPLIT_SCHUR_FACT_LOWER,
379: PC_FIELDSPLIT_SCHUR_FACT_UPPER,
380: PC_FIELDSPLIT_SCHUR_FACT_FULL
381: } PCFieldSplitSchurFactType;
382: PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
384: PETSC_EXTERN PETSC_DEPRECATED("Use PCFieldSplitSetSchurPre") PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
385: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC,PCFieldSplitSchurPreType,Mat);
386: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC,PCFieldSplitSchurPreType*,Mat*);
387: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
388: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
389: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC,Mat *S);
390: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC,Mat *S);
392: PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
393: PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
395: PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
396: PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);
398: PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);
400: PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
401: PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);
403: PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
404: PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);
406: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
407: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
408: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
410: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
411: PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
412: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
413: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
414: /*E
415: PCPARMSGlobalType - Determines the global preconditioner method in PARMS
417: Level: intermediate
419: .seealso: PCPARMSSetGlobal()
420: E*/
421: typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
422: PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
423: /*E
424: PCPARMSLocalType - Determines the local preconditioner method in PARMS
426: Level: intermediate
428: .seealso: PCPARMSSetLocal()
429: E*/
430: typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
431: PETSC_EXTERN const char *const PCPARMSLocalTypes[];
433: PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type);
434: PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type);
435: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits);
436: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart);
437: PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym);
438: PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2);
440: /*E
441: PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method
443: Level: intermediate
445: .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseProl()
446: E*/
447: typedef const char *PCGAMGType;
448: #define PCGAMGAGG "agg"
449: #define PCGAMGGEO "geo"
450: #define PCGAMGCLASSICAL "classical"
451: PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
452: PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
453: PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
454: PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
455: PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
456: PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
457: PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
458: PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType );
459: PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n);
460: PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n);
461: PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool);
462: PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool);
463: PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void);
464: PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void);
466: typedef const char *PCGAMGClassicalType;
467: #define PCGAMGCLASSICALDIRECT "direct"
468: #define PCGAMGCLASSICALSTANDARD "standard"
469: PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC,PCGAMGClassicalType);
471: #if defined(PETSC_HAVE_PCBDDC)
472: PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisLocalMat(PC,Mat);
473: PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC,IS);
474: PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
475: PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC,PetscInt);
476: PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace);
477: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
478: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC,IS);
479: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
480: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC,IS*);
481: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
482: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC,IS);
483: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
484: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC,IS*);
485: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
486: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC,PetscInt,IS[]);
487: PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode);
488: PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*);
489: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec);
490: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec);
491: #endif
493: PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool);
494: PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
495: PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec);
497: /*E
498: PCMGType - Determines the type of multigrid method that is run.
500: Level: beginner
502: Values:
503: + PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles()
504: . PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
505: smoothed before updating the residual. This only uses the
506: down smoother, in the preconditioner the upper smoother is ignored
507: . PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
508: that is starts on the coarsest grid, performs a cycle, interpolates
509: to the next, performs a cycle etc. This is much like the F-cycle presented in "Multigrid" by Trottenberg, Oosterlee, Schuller page 49, but that
510: algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
511: calls the V-cycle only on the coarser level and has a post-smoother instead.
512: - PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
513: from a finer
515: .seealso: PCMGSetType()
517: E*/
518: typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
519: PETSC_EXTERN const char *const PCMGTypes[];
520: #define PC_MG_CASCADE PC_MG_KASKADE;
522: /*E
523: PCMGCycleType - Use V-cycle or W-cycle
525: Level: beginner
527: Values:
528: + PC_MG_V_CYCLE
529: - PC_MG_W_CYCLE
531: .seealso: PCMGSetCycleType()
533: E*/
534: typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
535: PETSC_EXTERN const char *const PCMGCycleTypes[];
537: PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType);
538: PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*);
539: PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*);
541: PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt);
542: PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt);
543: PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType);
544: PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType);
545: PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt);
546: PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt);
547: PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool);
548: PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*);
551: PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec);
552: PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec);
553: PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec);
555: PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat);
556: PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*);
557: PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat);
558: PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*);
559: PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec);
560: PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*);
561: PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
562: PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat,Vec,Vec,Vec);
564: /*E
565: PCExoticType - Face based or wirebasket based coarse grid space
567: Level: beginner
569: .seealso: PCExoticSetType(), PCEXOTIC
570: E*/
571: typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
572: PETSC_EXTERN const char *const PCExoticTypes[];
573: PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);
575: #endif /* __PETSCPC_H */