Actual source code: petscpc.h
1: /* $Id: petscpc.h,v 1.122 2001/08/21 21:03:12 bsmith Exp $ */
3: /*
4: Preconditioner module.
5: */
8: #include petscmat.h
10: /*
11: PCList contains the list of preconditioners currently registered
12: These are added with the PCRegisterDynamic() macro
13: */
14: extern PetscFList PCList;
15: typedef char *PCType;
17: /*S
18: PC - Abstract PETSc object that manages all preconditioners
20: Level: beginner
22: Concepts: preconditioners
24: .seealso: PCCreate(), PCSetType(), PCType
25: S*/
26: typedef struct _p_PC* PC;
28: /*E
29: PCType - String with the name of a PETSc preconditioner method or the creation function
30: with an optional dynamic library name, for example
31: http://www.mcs.anl.gov/petsc/lib.a:mypccreate()
33: Level: beginner
35: .seealso: PCSetType(), PC
36: E*/
37: #define PCNONE "none"
38: #define PCJACOBI "jacobi"
39: #define PCSOR "sor"
40: #define PCLU "lu"
41: #define PCSHELL "shell"
42: #define PCBJACOBI "bjacobi"
43: #define PCMG "mg"
44: #define PCEISENSTAT "eisenstat"
45: #define PCILU "ilu"
46: #define PCICC "icc"
47: #define PCASM "asm"
48: #define PCSLES "sles"
49: #define PCCOMPOSITE "composite"
50: #define PCREDUNDANT "redundant"
51: #define PCSPAI "spai"
52: #define PCMILU "milu"
53: #define PCNN "nn"
54: #define PCCHOLESKY "cholesky"
55: #define PCRAMG "ramg"
56: #define PCPBJACOBI "pbjacobi"
57: #define PCMULTILEVEL "multilevel"
58: #define PCSCHUR "schur"
59: #define PCESI "esi"
60: #define PCPETSCESI "petscesi"
62: /* Logging support */
63: extern int PC_COOKIE;
64: extern int PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
65: extern int PC_ApplySymmetricRight, PC_ModifySubMatrices;
67: /*E
68: PCSide - If the preconditioner is to be applied to the left, right
69: or symmetrically around the operator.
71: Level: beginner
73: .seealso:
74: E*/
75: typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
77: EXTERN int PCCreate(MPI_Comm,PC*);
78: EXTERN int PCSetType(PC,PCType);
79: EXTERN int PCSetUp(PC);
80: EXTERN int PCSetUpOnBlocks(PC);
81: EXTERN int PCApply(PC,Vec,Vec);
82: EXTERN int PCApplySymmetricLeft(PC,Vec,Vec);
83: EXTERN int PCApplySymmetricRight(PC,Vec,Vec);
84: EXTERN int PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
85: EXTERN int PCApplyTranspose(PC,Vec,Vec);
86: EXTERN int PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
87: EXTERN int PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,int);
88: EXTERN int PCApplyRichardsonExists(PC,PetscTruth*);
90: EXTERN int PCRegisterDestroy(void);
91: EXTERN int PCRegisterAll(char*);
92: extern PetscTruth PCRegisterAllCalled;
94: EXTERN int PCRegister(char*,char*,char*,int(*)(PC));
95: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
96: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
97: #else
98: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
99: #endif
101: EXTERN int PCDestroy(PC);
102: EXTERN int PCSetFromOptions(PC);
103: EXTERN int PCGetType(PC,PCType*);
105: EXTERN int PCGetFactoredMatrix(PC,Mat*);
106: EXTERN int PCSetModifySubMatrices(PC,int(*)(PC,int,IS*,IS*,Mat*,void*),void*);
107: EXTERN int PCModifySubMatrices(PC,int,IS*,IS*,Mat*,void*);
109: EXTERN int PCSetOperators(PC,Mat,Mat,MatStructure);
110: EXTERN int PCGetOperators(PC,Mat*,Mat*,MatStructure*);
112: EXTERN int PCSetVector(PC,Vec);
113: EXTERN int PCGetVector(PC,Vec*);
114: EXTERN int PCView(PC,PetscViewer);
116: EXTERN int PCSetOptionsPrefix(PC,char*);
117: EXTERN int PCAppendOptionsPrefix(PC,char*);
118: EXTERN int PCGetOptionsPrefix(PC,char**);
120: EXTERN int PCNullSpaceAttach(PC,MatNullSpace);
122: EXTERN int PCComputeExplicitOperator(PC,Mat*);
124: /*
125: These are used to provide extra scaling of preconditioned
126: operator for time-stepping schemes like in PVODE
127: */
128: EXTERN int PCDiagonalScale(PC,PetscTruth*);
129: EXTERN int PCDiagonalScaleLeft(PC,Vec,Vec);
130: EXTERN int PCDiagonalScaleRight(PC,Vec,Vec);
131: EXTERN int PCDiagonalScaleSet(PC,Vec);
133: /* ------------- options specific to particular preconditioners --------- */
135: EXTERN int PCJacobiSetUseRowMax(PC);
136: EXTERN int PCSORSetSymmetric(PC,MatSORType);
137: EXTERN int PCSORSetOmega(PC,PetscReal);
138: EXTERN int PCSORSetIterations(PC,int,int);
140: EXTERN int PCEisenstatSetOmega(PC,PetscReal);
141: EXTERN int PCEisenstatNoDiagonalScaling(PC);
143: #define USE_PRECONDITIONER_MATRIX 0
144: #define USE_TRUE_MATRIX 1
145: EXTERN int PCBJacobiSetUseTrueLocal(PC);
146: EXTERN int PCBJacobiSetTotalBlocks(PC,int,int*);
147: EXTERN int PCBJacobiSetLocalBlocks(PC,int,int*);
149: EXTERN int PCSLESSetUseTrue(PC);
151: EXTERN int PCShellSetApply(PC,int (*)(void*,Vec,Vec),void*);
152: EXTERN int PCShellSetSetUp(PC,int (*)(void*));
153: EXTERN int PCShellSetApplyRichardson(PC,int (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,int),void*);
154: EXTERN int PCShellSetView(PC,int (*)(void*,PetscViewer));
155: EXTERN int PCShellSetName(PC,char*);
156: EXTERN int PCShellGetName(PC,char**);
158: EXTERN int PCLUSetMatOrdering(PC,MatOrderingType);
159: EXTERN int PCLUSetReuseOrdering(PC,PetscTruth);
160: EXTERN int PCLUSetReuseFill(PC,PetscTruth);
161: EXTERN int PCLUSetUseInPlace(PC);
162: EXTERN int PCLUSetFill(PC,PetscReal);
163: EXTERN int PCLUSetDamping(PC,PetscReal);
164: EXTERN int PCLUSetPivoting(PC,PetscReal);
166: EXTERN int PCCholeskySetMatOrdering(PC,MatOrderingType);
167: EXTERN int PCCholeskySetReuseOrdering(PC,PetscTruth);
168: EXTERN int PCCholeskySetReuseFill(PC,PetscTruth);
169: EXTERN int PCCholeskySetUseInPlace(PC);
170: EXTERN int PCCholeskySetFill(PC,PetscReal);
171: EXTERN int PCCholeskySetDamping(PC,PetscReal);
173: EXTERN int PCILUSetMatOrdering(PC,MatOrderingType);
174: EXTERN int PCILUSetUseInPlace(PC);
175: EXTERN int PCILUSetFill(PC,PetscReal);
176: EXTERN int PCILUSetLevels(PC,int);
177: EXTERN int PCILUSetReuseOrdering(PC,PetscTruth);
178: EXTERN int PCILUSetUseDropTolerance(PC,PetscReal,PetscReal,int);
179: EXTERN int PCILUDTSetReuseFill(PC,PetscTruth);
180: EXTERN int PCILUSetAllowDiagonalFill(PC);
181: EXTERN int PCILUSetDamping(PC,PetscReal);
182: EXTERN int PCILUSetSinglePrecisionSolves(PC,PetscTruth);
184: EXTERN int PCICCSetMatOrdering(PC,MatOrderingType);
185: EXTERN int PCICCSetFill(PC,PetscReal);
186: EXTERN int PCICCSetLevels(PC,int);
188: EXTERN int PCASMSetLocalSubdomains(PC,int,IS *);
189: EXTERN int PCASMSetTotalSubdomains(PC,int,IS *);
190: EXTERN int PCASMSetOverlap(PC,int);
191: /*E
192: PCASMType - Type of additive Schwarz method to use
194: $ PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
195: $ and computed values in ghost regions are added together. Classical
196: $ standard additive Schwarz
197: $ PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
198: $ region are discarded. Default
199: $ PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
200: $ region are added back in
201: $ PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
202: $ not very good.
204: Level: beginner
206: .seealso: PCASMSetType()
207: E*/
208: typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
209: EXTERN int PCASMSetType(PC,PCASMType);
210: EXTERN int PCASMCreateSubdomains2D(int,int,int,int,int,int,int *,IS **);
211: EXTERN int PCASMSetUseInPlace(PC);
212: EXTERN int PCASMGetLocalSubdomains(PC,int*,IS**);
213: EXTERN int PCASMGetLocalSubmatrices(PC,int*,Mat**);
215: /*E
216: PCCompositeType - Determines how two or more preconditioner are composed
218: $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
219: $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
220: $ computed after the previous preconditioner application
221: $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
222: $ where first preconditioner is built from alpha I + S and second from
223: $ alpha I + R
225: Level: beginner
227: .seealso: PCCompositeSetType()
228: E*/
229: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType;
230: EXTERN int PCCompositeSetUseTrue(PC);
231: EXTERN int PCCompositeSetType(PC,PCCompositeType);
232: EXTERN int PCCompositeAddPC(PC,PCType);
233: EXTERN int PCCompositeGetPC(PC pc,int n,PC *);
234: EXTERN int PCCompositeSpecialSetAlpha(PC,PetscScalar);
236: EXTERN int PCRedundantSetScatter(PC,VecScatter,VecScatter);
237: EXTERN int PCRedundantGetOperators(PC,Mat*,Mat*);
238: EXTERN int PCRedundantGetPC(PC,PC*);
239: EXTERN int MatGetOrderingList(PetscFList *list);
241: EXTERN int PCMultiLevelSetFields(PC, int, int);
242: EXTERN int PCMultiLevelSetNonlinearIterate(PC, Vec);
243: EXTERN int PCMultiLevelSetGradientOperator(PC, int, int, PetscScalar);
244: EXTERN int PCMultiLevelApplyGradient(PC, Vec, Vec);
245: EXTERN int PCMultiLevelApplyGradientTrans(PC, Vec, Vec);
246: EXTERN int PCMultiLevelBuildSolution(PC, Vec);
247: EXTERN int PCMultiLevelGetMultiplier(PC, Vec, Vec);
249: EXTERN int PCSchurSetGradientOperator(PC, int, int);
250: EXTERN int PCSchurGetIterationNumber(PC, int *, int *);
252: #endif /* __PETSCPC_H */