Actual source code: petscpc.h
petsc-3.4.2 2013-07-02
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 PCHMPI "hmpi"
68: #define PCSUPPORTGRAPH "supportgraph"
69: #define PCASA "asa"
70: #define PCCP "cp"
71: #define PCBFBT "bfbt"
72: #define PCLSC "lsc"
73: #define PCPYTHON "python"
74: #define PCPFMG "pfmg"
75: #define PCSYSPFMG "syspfmg"
76: #define PCREDISTRIBUTE "redistribute"
77: #define PCSVD "svd"
78: #define PCGAMG "gamg"
79: #define PCSACUSP "sacusp" /* these four run on NVIDIA GPUs using CUSP */
80: #define PCSACUSPPOLY "sacusppoly"
81: #define PCBICGSTABCUSP "bicgstabcusp"
82: #define PCAINVCUSP "ainvcusp"
83: #define PCBDDC "bddc"
85: /* Logging support */
86: PETSC_EXTERN PetscClassId PC_CLASSID;
88: /*E
89: PCSide - If the preconditioner is to be applied to the left, right
90: or symmetrically around the operator.
92: Level: beginner
94: .seealso:
95: E*/
96: typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
97: #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
98: PETSC_EXTERN const char *const *const PCSides;
100: PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
101: PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType);
102: PETSC_EXTERN PetscErrorCode PCSetUp(PC);
103: PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
104: PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
105: PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
106: PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
107: PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
108: PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
109: PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *);
110: PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
112: #define PC_FILE_CLASSID 1211222
114: /*E
115: PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
117: Level: advanced
119: Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
121: .seealso: PCApplyRichardson()
122: E*/
123: typedef enum {
124: PCRICHARDSON_CONVERGED_RTOL = 2,
125: PCRICHARDSON_CONVERGED_ATOL = 3,
126: PCRICHARDSON_CONVERGED_ITS = 4,
127: PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason;
129: PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
130: PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *);
131: PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool );
132: PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool);
133: PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*);
135: PETSC_EXTERN PetscErrorCode PCRegisterAll(void);
136: PETSC_EXTERN PetscBool PCRegisterAllCalled;
138: PETSC_EXTERN PetscErrorCode PCRegister(const char[],PetscErrorCode(*)(PC));
140: PETSC_EXTERN PetscErrorCode PCReset(PC);
141: PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
142: PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);
143: PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*);
145: PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
146: PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
147: PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
149: PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
150: PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
151: PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
153: PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer);
154: PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer);
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 PCFactorSetLevels(PC,PetscInt);
219: PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
221: PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
222: PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
223: PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
224: PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool);
225: PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*);
226: PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool );
228: /*E
229: PCASMType - Type of additive Schwarz method to use
231: $ PC_ASM_BASIC - Symmetric version where residuals from the ghost points are used
232: $ and computed values in ghost regions are added together.
233: $ Classical standard additive Schwarz.
234: $ PC_ASM_RESTRICT - Residuals from ghost points are used but computed values in ghost
235: $ region are discarded.
236: $ Default.
237: $ PC_ASM_INTERPOLATE - Residuals from ghost points are not used, computed values in ghost
238: $ region are added back in.
239: $ PC_ASM_NONE - Residuals from ghost points are not used, computed ghost values are
240: $ discarded.
241: $ Not very good.
243: Level: beginner
245: .seealso: PCASMSetType()
246: E*/
247: typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
248: PETSC_EXTERN const char *const PCASMTypes[];
250: PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
251: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
252: PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
253: PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
254: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
255: PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
257: /*E
258: PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
260: Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational
261: domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute
262: a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
263: (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
265: $ PC_GASM_BASIC - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
266: $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections
267: $ from neighboring subdomains.
268: $ Classical standard additive Schwarz.
269: $ PC_GASM_RESTRICT - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
270: $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result,
271: $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
272: $ assumption).
273: $ Default.
274: $ PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
275: $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
276: $ from neighboring subdomains.
277: $
278: $ PC_GASM_NONE - Residuals and corrections are zeroed out outside the local subdomains.
279: $ Not very good.
281: Level: beginner
283: .seealso: PCGASMSetType()
284: E*/
285: typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
286: PETSC_EXTERN const char *const PCGASMTypes[];
288: PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
289: PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool);
290: PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
291: PETSC_EXTERN PetscErrorCode PCGASMSetDMSubdomains(PC,PetscBool);
292: PETSC_EXTERN PetscErrorCode PCGASMGetDMSubdomains(PC,PetscBool*);
293: PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );
295: PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
296: PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]);
297: PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
298: PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
299: PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
300: PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);
302: /*E
303: PCCompositeType - Determines how two or more preconditioner are composed
305: $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
306: $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
307: $ computed after the previous preconditioner application
308: $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
309: $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
310: $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
311: $ where first preconditioner is built from alpha I + S and second from
312: $ alpha I + R
314: Level: beginner
316: .seealso: PCCompositeSetType()
317: E*/
318: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
319: PETSC_EXTERN const char *const PCCompositeTypes[];
321: PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
322: PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
323: PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
324: PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
326: PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
327: PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
328: PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
330: PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
331: PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
332: PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
333: PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
334: PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
335: PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
336: PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
337: PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
339: PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
340: PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
341: PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
342: PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
344: PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
345: PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
346: PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
347: PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
348: PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
349: PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
350: PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool);
351: PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*);
353: /*E
354: PCFieldSplitSchurPreType - Determines how to precondition Schur complement
356: Level: intermediate
358: .seealso: PCFieldSplitSchurPrecondition()
359: E*/
360: typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_A11,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
361: PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
363: /*E
364: PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
366: Level: intermediate
368: .seealso: PCFieldSplitSetSchurFactType()
369: E*/
370: typedef enum {
371: PC_FIELDSPLIT_SCHUR_FACT_DIAG,
372: PC_FIELDSPLIT_SCHUR_FACT_LOWER,
373: PC_FIELDSPLIT_SCHUR_FACT_UPPER,
374: PC_FIELDSPLIT_SCHUR_FACT_FULL
375: } PCFieldSplitSchurFactType;
376: PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
378: PETSC_EXTERN PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
379: PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
380: PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
382: PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
383: PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
385: PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
386: PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);
388: PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);
390: PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
391: PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);
393: PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
394: PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);
396: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
397: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
398: PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
400: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
401: PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
402: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
403: PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
404: /*E
405: PCPARMSGlobalType - Determines the global preconditioner method in PARMS
407: Level: intermediate
409: .seealso: PCPARMSSetGlobal()
410: E*/
411: typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
412: PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
413: /*E
414: PCPARMSLocalType - Determines the local preconditioner method in PARMS
416: Level: intermediate
418: .seealso: PCPARMSSetLocal()
419: E*/
420: typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
421: PETSC_EXTERN const char *const PCPARMSLocalTypes[];
423: PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type);
424: PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type);
425: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits);
426: PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart);
427: PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym);
428: PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2);
430: typedef const char *PCGAMGType;
431: #define PCGAMGAGG "agg"
432: #define PCGAMGGEO "geo"
433: PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
434: PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
435: PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
436: PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
437: PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
438: PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
439: PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
440: PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType );
441: PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n);
442: PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n);
443: PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool);
444: PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool);
445: PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void);
446: PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void);
448: #if defined(PETSC_HAVE_PCBDDC)
449: /* Enum defining how to treat the coarse problem */
450: typedef enum {SEQUENTIAL_BDDC,REPLICATED_BDDC,PARALLEL_BDDC,MULTILEVEL_BDDC} CoarseProblemType;
451: PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
452: PETSC_EXTERN PetscErrorCode PCBDDCSetMaxLevels(PC,PetscInt);
453: PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace);
454: PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
455: PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
456: PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
457: PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
458: PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseProblemType(PC,CoarseProblemType);
459: PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
460: PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode);
461: PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*);
462: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec);
463: PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec);
464: #endif
466: PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool);
467: PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
468: PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec);
470: /*E
471: PCMGType - Determines the type of multigrid method that is run.
473: Level: beginner
475: Values:
476: + PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles()
477: . PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
478: smoothed before updating the residual. This only uses the
479: down smoother, in the preconditioner the upper smoother is ignored
480: . PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
481: that is starts on the coarsest grid, performs a cycle, interpolates
482: 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
483: algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
484: calls the V-cycle only on the coarser level and has a post-smoother instead.
485: - PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
486: from a finer
488: .seealso: PCMGSetType()
490: E*/
491: typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
492: PETSC_EXTERN const char *const PCMGTypes[];
493: #define PC_MG_CASCADE PC_MG_KASKADE;
495: /*E
496: PCMGCycleType - Use V-cycle or W-cycle
498: Level: beginner
500: Values:
501: + PC_MG_V_CYCLE
502: - PC_MG_W_CYCLE
504: .seealso: PCMGSetCycleType()
506: E*/
507: typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
508: PETSC_EXTERN const char *const PCMGCycleTypes[];
510: PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType);
511: PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*);
512: PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*);
514: PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt);
515: PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt);
516: PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType);
517: PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType);
518: PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt);
519: PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt);
520: PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool);
521: PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*);
524: PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec);
525: PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec);
526: PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec);
528: PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat);
529: PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*);
530: PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat);
531: PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*);
532: PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec);
533: PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*);
534: PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
536: /*E
537: PCExoticType - Face based or wirebasket based coarse grid space
539: Level: beginner
541: .seealso: PCExoticSetType(), PCEXOTIC
542: E*/
543: typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
544: PETSC_EXTERN const char *const PCExoticTypes[];
545: PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);
547: #endif /* __PETSCPC_H */