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