xref: /petsc/include/petscpc.h (revision 2cfcea2917e9827eeaf4fa9c2e8dd75b054ced7b)
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 PETSCKSP_DLLEXPORT  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 const 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 PCSAMG            "samg"
58 #define PCPBJACOBI        "pbjacobi"
59 #define PCMAT             "mat"
60 #define PCHYPRE           "hypre"
61 #define PCFIELDSPLIT      "fieldsplit"
62 #define PCTFS             "tfs"
63 #define PCML              "ml"
64 #define PCPROMETHEUS      "prometheus"
65 
66 /* Logging support */
67 extern PetscCookie PETSCKSP_DLLEXPORT PC_COOKIE;
68 
69 /*E
70     PCSide - If the preconditioner is to be applied to the left, right
71      or symmetrically around the operator.
72 
73    Level: beginner
74 
75 .seealso:
76 E*/
77 typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
78 extern const char *PCSides[];
79 
80 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCreate(MPI_Comm,PC*);
81 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetType(PC,PCType);
82 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetUp(PC);
83 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetUpOnBlocks(PC);
84 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApply(PC,Vec,Vec);
85 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricLeft(PC,Vec,Vec);
86 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricRight(PC,Vec,Vec);
87 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
88 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyTranspose(PC,Vec,Vec);
89 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCHasApplyTranspose(PC,PetscTruth*);
90 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
91 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt);
92 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardsonExists(PC,PetscTruth*);
93 
94 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterDestroy(void);
95 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterAll(const char[]);
96 extern PetscTruth PCRegisterAllCalled;
97 
98 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
99 
100 /*MC
101    PCRegisterDynamic - Adds a method to the preconditioner package.
102 
103    Synopsis:
104    PetscErrorCode PCRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(PC))
105 
106    Not collective
107 
108    Input Parameters:
109 +  name_solver - name of a new user-defined solver
110 .  path - path (either absolute or relative) the library containing this solver
111 .  name_create - name of routine to create method context
112 -  routine_create - routine to create method context
113 
114    Notes:
115    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
116 
117    If dynamic libraries are used, then the fourth input argument (routine_create)
118    is ignored.
119 
120    Sample usage:
121 .vb
122    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
123               "MySolverCreate",MySolverCreate);
124 .ve
125 
126    Then, your solver can be chosen with the procedural interface via
127 $     PCSetType(pc,"my_solver")
128    or at runtime via the option
129 $     -pc_type my_solver
130 
131    Level: advanced
132 
133    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
134            occuring in pathname will be replaced with appropriate values.
135          If your function is not being put into a shared library then use PCRegister() instead
136 
137 .keywords: PC, register
138 
139 .seealso: PCRegisterAll(), PCRegisterDestroy()
140 M*/
141 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
142 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
143 #else
144 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
145 #endif
146 
147 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDestroy(PC);
148 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetFromOptions(PC);
149 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetType(PC,PCType*);
150 
151 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetFactoredMatrix(PC,Mat*);
152 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
153 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
154 
155 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetOperators(PC,Mat,Mat,MatStructure);
156 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperators(PC,Mat*,Mat*,MatStructure*);
157 
158 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCView(PC,PetscViewer);
159 
160 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetOptionsPrefix(PC,const char[]);
161 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCAppendOptionsPrefix(PC,const char[]);
162 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOptionsPrefix(PC,const char*[]);
163 
164 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCComputeExplicitOperator(PC,Mat*);
165 
166 /*
167       These are used to provide extra scaling of preconditioned
168    operator for time-stepping schemes like in SUNDIALS
169 */
170 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScale(PC,PetscTruth*);
171 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleLeft(PC,Vec,Vec);
172 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleRight(PC,Vec,Vec);
173 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleSet(PC,Vec);
174 
175 /* ------------- options specific to particular preconditioners --------- */
176 
177 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowMax(PC);
178 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseAbs(PC);
179 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric(PC,MatSORType);
180 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega(PC,PetscReal);
181 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations(PC,PetscInt,PetscInt);
182 
183 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatSetOmega(PC,PetscReal);
184 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatNoDiagonalScaling(PC);
185 
186 #define USE_PRECONDITIONER_MATRIX 0
187 #define USE_TRUE_MATRIX           1
188 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetUseTrueLocal(PC);
189 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
190 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
191 
192 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPSetUseTrue(PC);
193 
194 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApply(PC,PetscErrorCode (*)(void*,Vec,Vec));
195 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyTranspose(PC,PetscErrorCode (*)(void*,Vec,Vec));
196 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetSetUp(PC,PetscErrorCode (*)(void*));
197 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyRichardson(PC,PetscErrorCode (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt));
198 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetView(PC,PetscErrorCode (*)(void*,PetscViewer));
199 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetDestroy(PC,PetscErrorCode (*)(void*));
200 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetContext(PC,void**);
201 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetContext(PC,void*);
202 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetName(PC,const char[]);
203 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetName(PC,char*[]);
204 
205 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot(PC,PetscReal);
206 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero(PC,PetscReal);
207 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd(PC,PetscTruth);
208 
209 
210 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC,PetscReal);
211 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC,PetscReal);
212 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
213 
214 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrdering(PC,MatOrderingType);
215 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC,PetscTruth);
216 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC,PetscTruth);
217 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC);
218 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC);
219 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC,PetscTruth);
220 
221 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC,PetscInt);
222 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC,PetscReal,PetscReal,PetscInt);
223 
224 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetLocalSubdomains(PC,PetscInt,IS[]);
225 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetTotalSubdomains(PC,PetscInt,IS[]);
226 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetOverlap(PC,PetscInt);
227 /*E
228     PCASMType - Type of additive Schwarz method to use
229 
230 $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
231 $                 and computed values in ghost regions are added together. Classical
232 $                 standard additive Schwarz
233 $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
234 $                    region are discarded. Default
235 $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
236 $                       region are added back in
237 $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
238 $                not very good.
239 
240    Level: beginner
241 
242 .seealso: PCASMSetType()
243 E*/
244 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
245 extern const char *PCASMTypes[];
246 
247 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetType(PC,PCASMType);
248 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt *,IS **);
249 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetUseInPlace(PC);
250 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubdomains(PC,PetscInt*,IS*[]);
251 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
252 
253 /*E
254     PCCompositeType - Determines how two or more preconditioner are composed
255 
256 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
257 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
258 $                                computed after the previous preconditioner application
259 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
260 $                         where first preconditioner is built from alpha I + S and second from
261 $                         alpha I + R
262 
263    Level: beginner
264 
265 .seealso: PCCompositeSetType()
266 E*/
267 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType;
268 extern const char *PCCompositeTypes[];
269 
270 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC);
271 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC,PCCompositeType);
272 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC,PCType);
273 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC pc,PetscInt n,PC *);
274 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC,PetscScalar);
275 
276 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter(PC,VecScatter,VecScatter);
277 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators(PC,Mat*,Mat*);
278 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC(PC,PC*);
279 
280 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetEpsilon(PC,double);
281 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetNBSteps(PC,PetscInt);
282 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMax(PC,PetscInt);
283 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMaxNew(PC,PetscInt);
284 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetBlockSize(PC,PetscInt);
285 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetCacheSize(PC,PetscInt);
286 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetVerbose(PC,PetscInt);
287 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetSp(PC,PetscInt);
288 
289 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC,const char[]);
290 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
291 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
292 
293 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC,PetscInt,PetscInt*);
294 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC,PCCompositeType);
295 
296 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetCoordinates(PC,PetscReal*);
297 
298 PETSC_EXTERN_CXX_END
299 #endif /* __PETSCPC_H */
300