xref: /petsc/include/petscpc.h (revision 8229bfc2e272e22dd52541101af2fe6fdaaa77c5)
1 /* $Id: petscpc.h,v 1.122 2001/08/21 21:03:12 bsmith Exp $ */
2 
3 /*
4       Preconditioner module.
5 */
6 #if !defined(__PETSCPC_H)
7 #define __PETSCPC_H
8 #include "petscmat.h"
9 
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;
16 
17 /*S
18      PC - Abstract PETSc object that manages all preconditioners
19 
20    Level: beginner
21 
22   Concepts: preconditioners
23 
24 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
25 S*/
26 typedef struct _p_PC* PC;
27 
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()
32 
33    Level: beginner
34 
35    Notes: Click on the links below to see details on a particular solver
36 
37 .seealso: PCSetType(), PC, PCCreate()
38 E*/
39 #define PCNONE      "none"
40 #define PCJACOBI    "jacobi"
41 #define PCSOR       "sor"
42 #define PCLU        "lu"
43 #define PCSHELL     "shell"
44 #define PCBJACOBI   "bjacobi"
45 #define PCMG        "mg"
46 #define PCEISENSTAT "eisenstat"
47 #define PCILU       "ilu"
48 #define PCICC       "icc"
49 #define PCASM       "asm"
50 #define PCSLES      "sles"
51 #define PCCOMPOSITE "composite"
52 #define PCREDUNDANT "redundant"
53 #define PCSPAI      "spai"
54 #define PCMILU      "milu"
55 #define PCNN        "nn"
56 #define PCCHOLESKY  "cholesky"
57 #define PCRAMG      "ramg"
58 #define PCSAMG      "samg"
59 #define PCPBJACOBI  "pbjacobi"
60 #define PCMULTILEVEL "multilevel"
61 #define PCSCHUR      "schur"
62 #define PCESI        "esi"
63 #define PCPETSCESI   "petscesi"
64 #define PCMAT        "mat"
65 #define PCHYPRE      "hypre"
66 
67 /* Logging support */
68 extern int PC_COOKIE;
69 extern int PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
70 extern int PC_ApplySymmetricRight, PC_ModifySubMatrices;
71 
72 /*E
73     PCSide - If the preconditioner is to be applied to the left, right
74      or symmetrically around the operator.
75 
76    Level: beginner
77 
78 .seealso:
79 E*/
80 typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
81 
82 EXTERN int PCCreate(MPI_Comm,PC*);
83 EXTERN int PCSetType(PC,PCType);
84 EXTERN int PCSetUp(PC);
85 EXTERN int PCSetUpOnBlocks(PC);
86 EXTERN int PCApply(PC,Vec,Vec);
87 EXTERN int PCApplySymmetricLeft(PC,Vec,Vec);
88 EXTERN int PCApplySymmetricRight(PC,Vec,Vec);
89 EXTERN int PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
90 EXTERN int PCApplyTranspose(PC,Vec,Vec);
91 EXTERN int PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
92 EXTERN int PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,int);
93 EXTERN int PCApplyRichardsonExists(PC,PetscTruth*);
94 
95 EXTERN int        PCRegisterDestroy(void);
96 EXTERN int        PCRegisterAll(char*);
97 extern PetscTruth PCRegisterAllCalled;
98 
99 EXTERN int PCRegister(char*,char*,char*,int(*)(PC));
100 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
101 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
102 #else
103 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
104 #endif
105 
106 EXTERN int PCDestroy(PC);
107 EXTERN int PCSetFromOptions(PC);
108 EXTERN int PCESISetFromOptions(PC);
109 EXTERN int PCGetType(PC,PCType*);
110 
111 EXTERN int PCGetFactoredMatrix(PC,Mat*);
112 EXTERN int PCSetModifySubMatrices(PC,int(*)(PC,int,IS*,IS*,Mat*,void*),void*);
113 EXTERN int PCModifySubMatrices(PC,int,IS*,IS*,Mat*,void*);
114 
115 EXTERN int PCSetOperators(PC,Mat,Mat,MatStructure);
116 EXTERN int PCGetOperators(PC,Mat*,Mat*,MatStructure*);
117 
118 EXTERN int PCSetVector(PC,Vec);
119 EXTERN int PCGetVector(PC,Vec*);
120 EXTERN int PCView(PC,PetscViewer);
121 
122 EXTERN int PCSetOptionsPrefix(PC,char*);
123 EXTERN int PCAppendOptionsPrefix(PC,char*);
124 EXTERN int PCGetOptionsPrefix(PC,char**);
125 
126 EXTERN int PCNullSpaceAttach(PC,MatNullSpace);
127 
128 EXTERN int PCComputeExplicitOperator(PC,Mat*);
129 
130 /*
131       These are used to provide extra scaling of preconditioned
132    operator for time-stepping schemes like in PVODE
133 */
134 EXTERN int PCDiagonalScale(PC,PetscTruth*);
135 EXTERN int PCDiagonalScaleLeft(PC,Vec,Vec);
136 EXTERN int PCDiagonalScaleRight(PC,Vec,Vec);
137 EXTERN int PCDiagonalScaleSet(PC,Vec);
138 
139 /* ------------- options specific to particular preconditioners --------- */
140 
141 EXTERN int PCJacobiSetUseRowMax(PC);
142 EXTERN int PCSORSetSymmetric(PC,MatSORType);
143 EXTERN int PCSORSetOmega(PC,PetscReal);
144 EXTERN int PCSORSetIterations(PC,int,int);
145 
146 EXTERN int PCEisenstatSetOmega(PC,PetscReal);
147 EXTERN int PCEisenstatNoDiagonalScaling(PC);
148 
149 #define USE_PRECONDITIONER_MATRIX 0
150 #define USE_TRUE_MATRIX           1
151 EXTERN int PCBJacobiSetUseTrueLocal(PC);
152 EXTERN int PCBJacobiSetTotalBlocks(PC,int,int*);
153 EXTERN int PCBJacobiSetLocalBlocks(PC,int,int*);
154 
155 EXTERN int PCSLESSetUseTrue(PC);
156 
157 EXTERN int PCShellSetApply(PC,int (*)(void*,Vec,Vec),void*);
158 EXTERN int PCShellSetApplyTranspose(PC,int (*)(void*,Vec,Vec));
159 EXTERN int PCShellSetSetUp(PC,int (*)(void*));
160 EXTERN int PCShellSetApplyRichardson(PC,int (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,int),void*);
161 EXTERN int PCShellSetView(PC,int (*)(void*,PetscViewer));
162 EXTERN int PCShellSetName(PC,char*);
163 EXTERN int PCShellGetName(PC,char**);
164 
165 EXTERN int PCLUSetMatOrdering(PC,MatOrderingType);
166 EXTERN int PCLUSetReuseOrdering(PC,PetscTruth);
167 EXTERN int PCLUSetReuseFill(PC,PetscTruth);
168 EXTERN int PCLUSetUseInPlace(PC);
169 EXTERN int PCLUSetFill(PC,PetscReal);
170 EXTERN int PCLUSetDamping(PC,PetscReal);
171 EXTERN int PCLUSetPivoting(PC,PetscReal);
172 EXTERN int PCLUSetPivotInBlocks(PC,PetscTruth);
173 EXTERN int PCLUSetZeroPivot(PC,PetscReal);
174 
175 EXTERN int PCCholeskySetMatOrdering(PC,MatOrderingType);
176 EXTERN int PCCholeskySetReuseOrdering(PC,PetscTruth);
177 EXTERN int PCCholeskySetReuseFill(PC,PetscTruth);
178 EXTERN int PCCholeskySetUseInPlace(PC);
179 EXTERN int PCCholeskySetFill(PC,PetscReal);
180 EXTERN int PCCholeskySetDamping(PC,PetscReal);
181 EXTERN int PCCholeskySetPivotInBlocks(PC,PetscTruth);
182 
183 EXTERN int PCILUSetMatOrdering(PC,MatOrderingType);
184 EXTERN int PCILUSetUseInPlace(PC);
185 EXTERN int PCILUSetFill(PC,PetscReal);
186 EXTERN int PCILUSetLevels(PC,int);
187 EXTERN int PCILUSetReuseOrdering(PC,PetscTruth);
188 EXTERN int PCILUSetUseDropTolerance(PC,PetscReal,PetscReal,int);
189 EXTERN int PCILUDTSetReuseFill(PC,PetscTruth);
190 EXTERN int PCILUSetAllowDiagonalFill(PC);
191 EXTERN int PCILUSetDamping(PC,PetscReal);
192 EXTERN int PCILUSetPivotInBlocks(PC,PetscTruth);
193 EXTERN int PCILUSetZeroPivot(PC,PetscReal);
194 
195 EXTERN int PCICCSetMatOrdering(PC,MatOrderingType);
196 EXTERN int PCICCSetFill(PC,PetscReal);
197 EXTERN int PCICCSetLevels(PC,int);
198 EXTERN int PCICCSetPivotInBlocks(PC,PetscTruth);
199 
200 EXTERN int PCASMSetLocalSubdomains(PC,int,IS *);
201 EXTERN int PCASMSetTotalSubdomains(PC,int,IS *);
202 EXTERN int PCASMSetOverlap(PC,int);
203 /*E
204     PCASMType - Type of additive Schwarz method to use
205 
206 $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
207 $                 and computed values in ghost regions are added together. Classical
208 $                 standard additive Schwarz
209 $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
210 $                    region are discarded. Default
211 $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
212 $                       region are added back in
213 $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
214 $                not very good.
215 
216    Level: beginner
217 
218 .seealso: PCASMSetType()
219 E*/
220 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
221 EXTERN int PCASMSetType(PC,PCASMType);
222 EXTERN int PCASMCreateSubdomains2D(int,int,int,int,int,int,int *,IS **);
223 EXTERN int PCASMSetUseInPlace(PC);
224 EXTERN int PCASMGetLocalSubdomains(PC,int*,IS**);
225 EXTERN int PCASMGetLocalSubmatrices(PC,int*,Mat**);
226 
227 /*E
228     PCCompositeType - Determines how two or more preconditioner are composed
229 
230 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
231 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
232 $                                computed after the previous preconditioner application
233 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
234 $                         where first preconditioner is built from alpha I + S and second from
235 $                         alpha I + R
236 
237    Level: beginner
238 
239 .seealso: PCCompositeSetType()
240 E*/
241 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType;
242 EXTERN int PCCompositeSetUseTrue(PC);
243 EXTERN int PCCompositeSetType(PC,PCCompositeType);
244 EXTERN int PCCompositeAddPC(PC,PCType);
245 EXTERN int PCCompositeGetPC(PC pc,int n,PC *);
246 EXTERN int PCCompositeSpecialSetAlpha(PC,PetscScalar);
247 
248 EXTERN int PCRedundantSetScatter(PC,VecScatter,VecScatter);
249 EXTERN int PCRedundantGetOperators(PC,Mat*,Mat*);
250 EXTERN int PCRedundantGetPC(PC,PC*);
251 EXTERN int MatGetOrderingList(PetscFList *list);
252 
253 EXTERN int PCMultiLevelSetFields(PC, int, int);
254 EXTERN int PCMultiLevelSetNonlinearIterate(PC, Vec);
255 EXTERN int PCMultiLevelSetGradientOperator(PC, int, int, PetscScalar);
256 EXTERN int PCMultiLevelApplyGradient(PC, Vec, Vec);
257 EXTERN int PCMultiLevelApplyGradientTrans(PC, Vec, Vec);
258 EXTERN int PCMultiLevelBuildSolution(PC, Vec);
259 EXTERN int PCMultiLevelGetMultiplier(PC, Vec, Vec);
260 
261 EXTERN int PCSchurSetGradientOperator(PC, int, int);
262 EXTERN int PCSchurGetIterationNumber(PC, int *, int *);
263 
264 EXTERN int PCSPAISetEpsilon(PC,double);
265 EXTERN int PCSPAISetNBSteps(PC,int);
266 EXTERN int PCSPAISetMax(PC,int);
267 EXTERN int PCSPAISetMaxNew(PC,int);
268 EXTERN int PCSPAISetBlockSize(PC,int);
269 EXTERN int PCSPAISetCacheSize(PC,int);
270 EXTERN int PCSPAISetVerbose(PC,int);
271 EXTERN int PCSPAISetSp(PC,int);
272 
273 EXTERN int PCHYPRESetType(PC,char*);
274 EXTERN int PCBJacobiGetLocalBlocks(PC,int*,const int*[]);
275 EXTERN int PCBJacobiGetTotalBlocks(PC,int*,const int*[]);
276 
277 #endif /* __PETSCPC_H */
278