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