xref: /petsc/include/petscpc.h (revision 3b5a905bbcba14b9250de9c3c6bc067664eddbb4)
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
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 .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"
61 #define PCMAT        "mat"
62 #define PCHYPRE      "hypre"
63 
64 /* Logging support */
65 extern int PC_COOKIE;
66 extern int PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft;
67 extern int PC_ApplySymmetricRight, PC_ModifySubMatrices;
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 
79 EXTERN int PCCreate(MPI_Comm,PC*);
80 EXTERN int PCSetType(PC,PCType);
81 EXTERN int PCSetUp(PC);
82 EXTERN int PCSetUpOnBlocks(PC);
83 EXTERN int PCApply(PC,Vec,Vec);
84 EXTERN int PCApplySymmetricLeft(PC,Vec,Vec);
85 EXTERN int PCApplySymmetricRight(PC,Vec,Vec);
86 EXTERN int PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
87 EXTERN int PCApplyTranspose(PC,Vec,Vec);
88 EXTERN int PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
89 EXTERN int PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,int);
90 EXTERN int PCApplyRichardsonExists(PC,PetscTruth*);
91 
92 EXTERN int        PCRegisterDestroy(void);
93 EXTERN int        PCRegisterAll(char*);
94 extern PetscTruth PCRegisterAllCalled;
95 
96 EXTERN int PCRegister(char*,char*,char*,int(*)(PC));
97 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
98 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
99 #else
100 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
101 #endif
102 
103 EXTERN int PCDestroy(PC);
104 EXTERN int PCSetFromOptions(PC);
105 EXTERN int PCGetType(PC,PCType*);
106 
107 EXTERN int PCGetFactoredMatrix(PC,Mat*);
108 EXTERN int PCSetModifySubMatrices(PC,int(*)(PC,int,IS*,IS*,Mat*,void*),void*);
109 EXTERN int PCModifySubMatrices(PC,int,IS*,IS*,Mat*,void*);
110 
111 EXTERN int PCSetOperators(PC,Mat,Mat,MatStructure);
112 EXTERN int PCGetOperators(PC,Mat*,Mat*,MatStructure*);
113 
114 EXTERN int PCSetVector(PC,Vec);
115 EXTERN int PCGetVector(PC,Vec*);
116 EXTERN int PCView(PC,PetscViewer);
117 
118 EXTERN int PCSetOptionsPrefix(PC,char*);
119 EXTERN int PCAppendOptionsPrefix(PC,char*);
120 EXTERN int PCGetOptionsPrefix(PC,char**);
121 
122 EXTERN int PCNullSpaceAttach(PC,MatNullSpace);
123 
124 EXTERN int PCComputeExplicitOperator(PC,Mat*);
125 
126 /*
127       These are used to provide extra scaling of preconditioned
128    operator for time-stepping schemes like in PVODE
129 */
130 EXTERN int PCDiagonalScale(PC,PetscTruth*);
131 EXTERN int PCDiagonalScaleLeft(PC,Vec,Vec);
132 EXTERN int PCDiagonalScaleRight(PC,Vec,Vec);
133 EXTERN int PCDiagonalScaleSet(PC,Vec);
134 
135 /* ------------- options specific to particular preconditioners --------- */
136 
137 EXTERN int PCJacobiSetUseRowMax(PC);
138 EXTERN int PCSORSetSymmetric(PC,MatSORType);
139 EXTERN int PCSORSetOmega(PC,PetscReal);
140 EXTERN int PCSORSetIterations(PC,int,int);
141 
142 EXTERN int PCEisenstatSetOmega(PC,PetscReal);
143 EXTERN int PCEisenstatNoDiagonalScaling(PC);
144 
145 #define USE_PRECONDITIONER_MATRIX 0
146 #define USE_TRUE_MATRIX           1
147 EXTERN int PCBJacobiSetUseTrueLocal(PC);
148 EXTERN int PCBJacobiSetTotalBlocks(PC,int,int*);
149 EXTERN int PCBJacobiSetLocalBlocks(PC,int,int*);
150 
151 EXTERN int PCSLESSetUseTrue(PC);
152 
153 EXTERN int PCShellSetApply(PC,int (*)(void*,Vec,Vec),void*);
154 EXTERN int PCShellSetSetUp(PC,int (*)(void*));
155 EXTERN int PCShellSetApplyRichardson(PC,int (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,int),void*);
156 EXTERN int PCShellSetView(PC,int (*)(void*,PetscViewer));
157 EXTERN int PCShellSetName(PC,char*);
158 EXTERN int PCShellGetName(PC,char**);
159 
160 EXTERN int PCLUSetMatOrdering(PC,MatOrderingType);
161 EXTERN int PCLUSetReuseOrdering(PC,PetscTruth);
162 EXTERN int PCLUSetReuseFill(PC,PetscTruth);
163 EXTERN int PCLUSetUseInPlace(PC);
164 EXTERN int PCLUSetFill(PC,PetscReal);
165 EXTERN int PCLUSetDamping(PC,PetscReal);
166 EXTERN int PCLUSetPivoting(PC,PetscReal);
167 
168 EXTERN int PCCholeskySetMatOrdering(PC,MatOrderingType);
169 EXTERN int PCCholeskySetReuseOrdering(PC,PetscTruth);
170 EXTERN int PCCholeskySetReuseFill(PC,PetscTruth);
171 EXTERN int PCCholeskySetUseInPlace(PC);
172 EXTERN int PCCholeskySetFill(PC,PetscReal);
173 EXTERN int PCCholeskySetDamping(PC,PetscReal);
174 
175 EXTERN int PCILUSetMatOrdering(PC,MatOrderingType);
176 EXTERN int PCILUSetUseInPlace(PC);
177 EXTERN int PCILUSetFill(PC,PetscReal);
178 EXTERN int PCILUSetLevels(PC,int);
179 EXTERN int PCILUSetReuseOrdering(PC,PetscTruth);
180 EXTERN int PCILUSetUseDropTolerance(PC,PetscReal,PetscReal,int);
181 EXTERN int PCILUDTSetReuseFill(PC,PetscTruth);
182 EXTERN int PCILUSetAllowDiagonalFill(PC);
183 EXTERN int PCILUSetDamping(PC,PetscReal);
184 EXTERN int PCILUSetSinglePrecisionSolves(PC,PetscTruth);
185 
186 EXTERN int PCICCSetMatOrdering(PC,MatOrderingType);
187 EXTERN int PCICCSetFill(PC,PetscReal);
188 EXTERN int PCICCSetLevels(PC,int);
189 
190 EXTERN int PCASMSetLocalSubdomains(PC,int,IS *);
191 EXTERN int PCASMSetTotalSubdomains(PC,int,IS *);
192 EXTERN int PCASMSetOverlap(PC,int);
193 /*E
194     PCASMType - Type of additive Schwarz method to use
195 
196 $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
197 $                 and computed values in ghost regions are added together. Classical
198 $                 standard additive Schwarz
199 $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
200 $                    region are discarded. Default
201 $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
202 $                       region are added back in
203 $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
204 $                not very good.
205 
206    Level: beginner
207 
208 .seealso: PCASMSetType()
209 E*/
210 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
211 EXTERN int PCASMSetType(PC,PCASMType);
212 EXTERN int PCASMCreateSubdomains2D(int,int,int,int,int,int,int *,IS **);
213 EXTERN int PCASMSetUseInPlace(PC);
214 EXTERN int PCASMGetLocalSubdomains(PC,int*,IS**);
215 EXTERN int PCASMGetLocalSubmatrices(PC,int*,Mat**);
216 
217 /*E
218     PCCompositeType - Determines how two or more preconditioner are composed
219 
220 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
221 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
222 $                                computed after the previous preconditioner application
223 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
224 $                         where first preconditioner is built from alpha I + S and second from
225 $                         alpha I + R
226 
227    Level: beginner
228 
229 .seealso: PCCompositeSetType()
230 E*/
231 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType;
232 EXTERN int PCCompositeSetUseTrue(PC);
233 EXTERN int PCCompositeSetType(PC,PCCompositeType);
234 EXTERN int PCCompositeAddPC(PC,PCType);
235 EXTERN int PCCompositeGetPC(PC pc,int n,PC *);
236 EXTERN int PCCompositeSpecialSetAlpha(PC,PetscScalar);
237 
238 EXTERN int PCRedundantSetScatter(PC,VecScatter,VecScatter);
239 EXTERN int PCRedundantGetOperators(PC,Mat*,Mat*);
240 EXTERN int PCRedundantGetPC(PC,PC*);
241 EXTERN int MatGetOrderingList(PetscFList *list);
242 
243 EXTERN int PCMultiLevelSetFields(PC, int, int);
244 EXTERN int PCMultiLevelSetNonlinearIterate(PC, Vec);
245 EXTERN int PCMultiLevelSetGradientOperator(PC, int, int, PetscScalar);
246 EXTERN int PCMultiLevelApplyGradient(PC, Vec, Vec);
247 EXTERN int PCMultiLevelApplyGradientTrans(PC, Vec, Vec);
248 EXTERN int PCMultiLevelBuildSolution(PC, Vec);
249 EXTERN int PCMultiLevelGetMultiplier(PC, Vec, Vec);
250 
251 EXTERN int PCSchurSetGradientOperator(PC, int, int);
252 EXTERN int PCSchurGetIterationNumber(PC, int *, int *);
253 
254 #endif /* __PETSCPC_H */
255