xref: /petsc/include/petscpc.h (revision 7044107279799bb7fd2e4de09560a77cd45ddebe)
1 /*
2       Preconditioner module.
3 */
4 #if !defined(__PETSCPC_H)
5 #define __PETSCPC_H
6 #include "petscmat.h"
7 PETSC_EXTERN_CXX_BEGIN
8 
9 EXTERN PetscErrorCode PCInitializePackage(const char[]);
10 
11 /*
12     PCList contains the list of preconditioners currently registered
13    These are added with the PCRegisterDynamic() macro
14 */
15 extern PetscFList PCList;
16 #define PCType char*
17 
18 /*S
19      PC - Abstract PETSc object that manages all preconditioners
20 
21    Level: beginner
22 
23   Concepts: preconditioners
24 
25 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
26 S*/
27 typedef struct _p_PC* PC;
28 
29 /*E
30     PCType - String with the name of a PETSc preconditioner method or the creation function
31        with an optional dynamic library name, for example
32        http://www.mcs.anl.gov/petsc/lib.a:mypccreate()
33 
34    Level: beginner
35 
36    Notes: Click on the links below to see details on a particular solver
37 
38 .seealso: PCSetType(), PC, PCCreate()
39 E*/
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 PCKSP        "ksp"
52 #define PCCOMPOSITE  "composite"
53 #define PCREDUNDANT  "redundant"
54 #define PCSPAI       "spai"
55 #define PCNN         "nn"
56 #define PCCHOLESKY   "cholesky"
57 #define PCRAMG       "ramg"
58 #define PCSAMG       "samg"
59 #define PCPBJACOBI   "pbjacobi"
60 #define PCMAT        "mat"
61 #define PCHYPRE      "hypre"
62 
63 /* Logging support */
64 extern PetscCookie PC_COOKIE;
65 extern PetscEvent    PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
66 extern PetscEvent    PC_ApplySymmetricRight, PC_ModifySubMatrices;
67 
68 /*E
69     PCSide - If the preconditioner is to be applied to the left, right
70      or symmetrically around the operator.
71 
72    Level: beginner
73 
74 .seealso:
75 E*/
76 typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
77 
78 EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
79 EXTERN PetscErrorCode PCSetType(PC,const PCType);
80 EXTERN PetscErrorCode PCSetUp(PC);
81 EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
82 EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
83 EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
84 EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
85 EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
86 EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
87 EXTERN PetscErrorCode PCHasApplyTranspose(PC,PetscTruth*);
88 EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
89 EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,int);
90 EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscTruth*);
91 
92 EXTERN PetscErrorCode        PCRegisterDestroy(void);
93 EXTERN PetscErrorCode        PCRegisterAll(const char[]);
94 extern PetscTruth PCRegisterAllCalled;
95 
96 EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
97 
98 /*MC
99    PCRegisterDynamic - Adds a method to the preconditioner package.
100 
101    Synopsis:
102    PetscErrorCode PCRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(PC))
103 
104    Not collective
105 
106    Input Parameters:
107 +  name_solver - name of a new user-defined solver
108 .  path - path (either absolute or relative) the library containing this solver
109 .  name_create - name of routine to create method context
110 -  routine_create - routine to create method context
111 
112    Notes:
113    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
114 
115    If dynamic libraries are used, then the fourth input argument (routine_create)
116    is ignored.
117 
118    Sample usage:
119 .vb
120    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
121               "MySolverCreate",MySolverCreate);
122 .ve
123 
124    Then, your solver can be chosen with the procedural interface via
125 $     PCSetType(pc,"my_solver")
126    or at runtime via the option
127 $     -pc_type my_solver
128 
129    Level: advanced
130 
131    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT}, or ${any environmental variable}
132            occuring in pathname will be replaced with appropriate values.
133          If your function is not being put into a shared library then use PCRegister() instead
134 
135 .keywords: PC, register
136 
137 .seealso: PCRegisterAll(), PCRegisterDestroy()
138 M*/
139 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
140 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
141 #else
142 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
143 #endif
144 
145 EXTERN PetscErrorCode PCDestroy(PC);
146 EXTERN PetscErrorCode PCSetFromOptions(PC);
147 EXTERN PetscErrorCode PCGetType(PC,PCType*);
148 
149 EXTERN PetscErrorCode PCGetFactoredMatrix(PC,Mat*);
150 EXTERN PetscErrorCode PCSetModifySubMatrices(PC,int(*)(PC,int,const IS[],const IS[],Mat[],void*),void*);
151 EXTERN PetscErrorCode PCModifySubMatrices(PC,int,const IS[],const IS[],Mat[],void*);
152 
153 EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
154 EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
155 
156 EXTERN PetscErrorCode PCView(PC,PetscViewer);
157 
158 EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
159 EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
160 EXTERN PetscErrorCode PCGetOptionsPrefix(PC,char*[]);
161 
162 EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
163 
164 /*
165       These are used to provide extra scaling of preconditioned
166    operator for time-stepping schemes like in PVODE
167 */
168 EXTERN PetscErrorCode PCDiagonalScale(PC,PetscTruth*);
169 EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
170 EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
171 EXTERN PetscErrorCode PCDiagonalScaleSet(PC,Vec);
172 
173 /* ------------- options specific to particular preconditioners --------- */
174 
175 EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
176 EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
177 EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
178 EXTERN PetscErrorCode PCSORSetIterations(PC,int,int);
179 
180 EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
181 EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
182 
183 #define USE_PRECONDITIONER_MATRIX 0
184 #define USE_TRUE_MATRIX           1
185 EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC);
186 EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,int,const int[]);
187 EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,int,const int[]);
188 
189 EXTERN PetscErrorCode PCKSPSetUseTrue(PC);
190 
191 EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(void*,Vec,Vec),void*);
192 EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(void*,Vec,Vec));
193 EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(void*));
194 EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,int),void*);
195 EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(void*,PetscViewer));
196 EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
197 EXTERN PetscErrorCode PCShellGetName(PC,char*[]);
198 
199 EXTERN PetscErrorCode PCLUSetMatOrdering(PC,MatOrderingType);
200 EXTERN PetscErrorCode PCLUSetReuseOrdering(PC,PetscTruth);
201 EXTERN PetscErrorCode PCLUSetReuseFill(PC,PetscTruth);
202 EXTERN PetscErrorCode PCLUSetUseInPlace(PC);
203 EXTERN PetscErrorCode PCLUSetFill(PC,PetscReal);
204 EXTERN PetscErrorCode PCLUSetDamping(PC,PetscReal);
205 EXTERN PetscErrorCode PCLUSetShift(PC,PetscTruth);
206 EXTERN PetscErrorCode PCLUSetPivoting(PC,PetscReal);
207 EXTERN PetscErrorCode PCLUSetPivotInBlocks(PC,PetscTruth);
208 EXTERN PetscErrorCode PCLUSetZeroPivot(PC,PetscReal);
209 
210 EXTERN PetscErrorCode PCCholeskySetMatOrdering(PC,MatOrderingType);
211 EXTERN PetscErrorCode PCCholeskySetReuseOrdering(PC,PetscTruth);
212 EXTERN PetscErrorCode PCCholeskySetReuseFill(PC,PetscTruth);
213 EXTERN PetscErrorCode PCCholeskySetUseInPlace(PC);
214 EXTERN PetscErrorCode PCCholeskySetFill(PC,PetscReal);
215 EXTERN PetscErrorCode PCCholeskySetDamping(PC,PetscReal);
216 EXTERN PetscErrorCode PCCholeskySetShift(PC,PetscTruth);
217 EXTERN PetscErrorCode PCCholeskySetPivotInBlocks(PC,PetscTruth);
218 
219 EXTERN PetscErrorCode PCILUSetMatOrdering(PC,MatOrderingType);
220 EXTERN PetscErrorCode PCILUSetUseInPlace(PC);
221 EXTERN PetscErrorCode PCILUSetFill(PC,PetscReal);
222 EXTERN PetscErrorCode PCILUSetLevels(PC,PetscInt);
223 EXTERN PetscErrorCode PCILUSetReuseOrdering(PC,PetscTruth);
224 EXTERN PetscErrorCode PCILUSetUseDropTolerance(PC,PetscReal,PetscReal,PetscInt);
225 EXTERN PetscErrorCode PCILUDTSetReuseFill(PC,PetscTruth);
226 EXTERN PetscErrorCode PCILUSetAllowDiagonalFill(PC);
227 EXTERN PetscErrorCode PCILUSetDamping(PC,PetscReal);
228 EXTERN PetscErrorCode PCILUSetShift(PC,PetscTruth);
229 EXTERN PetscErrorCode PCILUSetPivotInBlocks(PC,PetscTruth);
230 EXTERN PetscErrorCode PCILUSetZeroPivot(PC,PetscReal);
231 
232 EXTERN PetscErrorCode PCICCSetMatOrdering(PC,MatOrderingType);
233 EXTERN PetscErrorCode PCICCSetFill(PC,PetscReal);
234 EXTERN PetscErrorCode PCICCSetLevels(PC,PetscInt);
235 EXTERN PetscErrorCode PCICCSetDamping(PC,PetscReal);
236 EXTERN PetscErrorCode PCICCSetShift(PC,PetscTruth);
237 EXTERN PetscErrorCode PCICCSetPivotInBlocks(PC,PetscTruth);
238 EXTERN PetscErrorCode PCICCSetZeroPivot(PC,PetscReal);
239 
240 EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[]);
241 EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[]);
242 EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
243 /*E
244     PCASMType - Type of additive Schwarz method to use
245 
246 $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
247 $                 and computed values in ghost regions are added together. Classical
248 $                 standard additive Schwarz
249 $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
250 $                    region are discarded. Default
251 $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
252 $                       region are added back in
253 $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
254 $                not very good.
255 
256    Level: beginner
257 
258 .seealso: PCASMSetType()
259 E*/
260 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
261 EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
262 EXTERN PetscErrorCode PCASMCreateSubdomains2D(int,int,int,int,int,int,int *,IS **);
263 EXTERN PetscErrorCode PCASMSetUseInPlace(PC);
264 EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,int*,IS*[]);
265 EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,int*,Mat*[]);
266 
267 /*E
268     PCCompositeType - Determines how two or more preconditioner are composed
269 
270 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
271 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
272 $                                computed after the previous preconditioner application
273 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
274 $                         where first preconditioner is built from alpha I + S and second from
275 $                         alpha I + R
276 
277    Level: beginner
278 
279 .seealso: PCCompositeSetType()
280 E*/
281 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType;
282 EXTERN PetscErrorCode PCCompositeSetUseTrue(PC);
283 EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
284 EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
285 EXTERN PetscErrorCode PCCompositeGetPC(PC pc,int n,PC *);
286 EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
287 
288 EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
289 EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
290 EXTERN PetscErrorCode PCRedundantGetPC(PC,PC*);
291 EXTERN PetscErrorCode MatGetOrderingList(PetscFList *list);
292 
293 EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
294 EXTERN PetscErrorCode PCSPAISetNBSteps(PC,int);
295 EXTERN PetscErrorCode PCSPAISetMax(PC,int);
296 EXTERN PetscErrorCode PCSPAISetMaxNew(PC,int);
297 EXTERN PetscErrorCode PCSPAISetBlockSize(PC,int);
298 EXTERN PetscErrorCode PCSPAISetCacheSize(PC,int);
299 EXTERN PetscErrorCode PCSPAISetVerbose(PC,int);
300 EXTERN PetscErrorCode PCSPAISetSp(PC,int);
301 
302 EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
303 EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,int*,const int*[]);
304 EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,int*,const int*[]);
305 
306 PETSC_EXTERN_CXX_END
307 #endif /* __PETSCPC_H */
308