xref: /petsc/include/petscpc.h (revision 1111b64dd641ff9b0189bccba1c2a1ee10f9b8d6)
1 /* $Id: petscpc.h,v 1.114 2001/01/29 19:25:52 bsmith Exp bsmith $ */
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 
18 #define PC_COOKIE     PETSC_COOKIE+9
19 /*S
20      PC - Abstract PETSc object that manages all preconditioners
21 
22    Level: beginner
23 
24   Concepts: preconditioners
25 
26 .seealso:  PCCreate(), PCSetType(), PCType
27 S*/
28 typedef struct _p_PC* PC;
29 
30 /*E
31     PCType - String with the name of a PETSc preconditioner method or the creation function
32        with an optional dynamic library name, for example
33        http://www.mcs.anl.gov/petsc/lib.a:mypccreate()
34 
35    Level: beginner
36 
37 .seealso: PCSetType(), PC
38 E*/
39 #define PCNONE      "none"
40 #define PCJACOBI    "jacobi"
41 #define PCSOR       "sor"
42 #define PCLU        "lu"
43 #define PCSHELL     "shell"
44 #define PCBJACOBI   "bjacobi"
45 #define PCMG        "mg"
46 #define PCEISENSTAT "eisenstat"
47 #define PCILU       "ilu"
48 #define PCICC       "icc"
49 #define PCASM       "asm"
50 #define PCSLES      "sles"
51 #define PCCOMPOSITE "composite"
52 #define PCREDUNDANT "redundant"
53 #define PCSPAI      "spai"
54 #define PCMILU      "milu"
55 #define PCNN        "nn"
56 #define PCCHOLESKY  "cholesky"
57 #define PCRAMG      "ramg"
58 
59 /*E
60     PCSide - If the preconditioner is to be applied to the left, right
61      or symmetrically around the operator.
62 
63    Level: beginner
64 
65 .seealso:
66 E*/
67 typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
68 
69 EXTERN int PCCreate(MPI_Comm,PC*);
70 EXTERN int PCSetType(PC,PCType);
71 EXTERN int PCSetUp(PC);
72 EXTERN int PCSetUpOnBlocks(PC);
73 EXTERN int PCApply(PC,Vec,Vec);
74 EXTERN int PCApplySymmetricLeft(PC,Vec,Vec);
75 EXTERN int PCApplySymmetricRight(PC,Vec,Vec);
76 EXTERN int PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
77 EXTERN int PCApplyTranspose(PC,Vec,Vec);
78 EXTERN int PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
79 EXTERN int PCApplyRichardson(PC,Vec,Vec,Vec,int);
80 EXTERN int PCApplyRichardsonExists(PC,PetscTruth*);
81 
82 EXTERN int        PCRegisterDestroy(void);
83 EXTERN int        PCRegisterAll(char*);
84 extern PetscTruth PCRegisterAllCalled;
85 
86 EXTERN int PCRegister(char*,char*,char*,int(*)(PC));
87 #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
88 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
89 #else
90 #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
91 #endif
92 
93 EXTERN int PCDestroy(PC);
94 EXTERN int PCSetFromOptions(PC);
95 EXTERN int PCSetTypeFromOptions(PC);
96 EXTERN int PCGetType(PC,PCType*);
97 
98 EXTERN int PCGetFactoredMatrix(PC,Mat*);
99 EXTERN int PCSetModifySubMatrices(PC,int(*)(PC,int,IS*,IS*,Mat*,void*),void*);
100 EXTERN int PCModifySubMatrices(PC,int,IS*,IS*,Mat*,void*);
101 
102 EXTERN int PCSetOperators(PC,Mat,Mat,MatStructure);
103 EXTERN int PCGetOperators(PC,Mat*,Mat*,MatStructure*);
104 
105 EXTERN int PCSetVector(PC,Vec);
106 EXTERN int PCGetVector(PC,Vec*);
107 EXTERN int PCView(PC,PetscViewer);
108 
109 EXTERN int PCSetOptionsPrefix(PC,char*);
110 EXTERN int PCAppendOptionsPrefix(PC,char*);
111 EXTERN int PCGetOptionsPrefix(PC,char**);
112 
113 EXTERN int PCNullSpaceAttach(PC,MatNullSpace);
114 
115 EXTERN int PCComputeExplicitOperator(PC,Mat*);
116 
117 /*
118       These are used to provide extra scaling of preconditioned
119    operator for time-stepping schemes like in PVODE
120 */
121 EXTERN int PCDiagonalScale(PC,PetscTruth*);
122 EXTERN int PCDiagonalScaleLeft(PC,Vec,Vec);
123 EXTERN int PCDiagonalScaleRight(PC,Vec,Vec);
124 EXTERN int PCDiagonalScaleSet(PC,Vec);
125 
126 /* ------------- options specific to particular preconditioners --------- */
127 
128 EXTERN int PCJacobiSetUseRowMax(PC);
129 EXTERN int PCSORSetSymmetric(PC,MatSORType);
130 EXTERN int PCSORSetOmega(PC,double);
131 EXTERN int PCSORSetIterations(PC,int);
132 
133 EXTERN int PCEisenstatSetOmega(PC,double);
134 EXTERN int PCEisenstatNoDiagonalScaling(PC);
135 
136 #define USE_PRECONDITIONER_MATRIX 0
137 #define USE_TRUE_MATRIX           1
138 EXTERN int PCBJacobiSetUseTrueLocal(PC);
139 EXTERN int PCBJacobiSetTotalBlocks(PC,int,int*);
140 EXTERN int PCBJacobiSetLocalBlocks(PC,int,int*);
141 
142 EXTERN int PCSLESSetUseTrue(PC);
143 EXTERN int PCCompositeSetUseTrue(PC);
144 
145 EXTERN int PCShellSetApply(PC,int (*)(void*,Vec,Vec),void*);
146 EXTERN int PCShellSetSetUp(PC,int (*)(void*));
147 EXTERN int PCShellSetApplyRichardson(PC,int (*)(void*,Vec,Vec,Vec,int),void*);
148 EXTERN int PCShellSetName(PC,char*);
149 EXTERN int PCShellGetName(PC,char**);
150 
151 EXTERN int PCLUSetMatOrdering(PC,MatOrderingType);
152 EXTERN int PCLUSetReuseOrdering(PC,PetscTruth);
153 EXTERN int PCLUSetReuseFill(PC,PetscTruth);
154 EXTERN int PCLUSetUseInPlace(PC);
155 EXTERN int PCLUSetFill(PC,double);
156 EXTERN int PCLUSetDamping(PC,double);
157 EXTERN int PCLUSetColumnPivoting(PC,PetscReal);
158 
159 EXTERN int PCCholeskySetMatOrdering(PC,MatOrderingType);
160 EXTERN int PCCholeskySetReuseOrdering(PC,PetscTruth);
161 EXTERN int PCCholeskySetReuseFill(PC,PetscTruth);
162 EXTERN int PCCholeskySetUseInPlace(PC);
163 EXTERN int PCCholeskySetFill(PC,double);
164 EXTERN int PCCholeskySetDamping(PC,double);
165 
166 EXTERN int PCILUSetMatOrdering(PC,MatOrderingType);
167 EXTERN int PCILUSetUseInPlace(PC);
168 EXTERN int PCILUSetFill(PC,double);
169 EXTERN int PCILUSetLevels(PC,int);
170 EXTERN int PCILUSetReuseOrdering(PC,PetscTruth);
171 EXTERN int PCILUSetUseDropTolerance(PC,PetscReal,PetscReal,int);
172 EXTERN int PCILUDTSetReuseFill(PC,PetscTruth);
173 EXTERN int PCILUSetAllowDiagonalFill(PC);
174 EXTERN int PCILUSetDamping(PC,double);
175 
176 EXTERN int PCICCSetMatOrdering(PC,MatOrderingType);
177 EXTERN int PCICCSetFill(PC,double);
178 EXTERN int PCICCSetLevels(PC,int);
179 
180 EXTERN int PCASMSetLocalSubdomains(PC,int,IS *);
181 EXTERN int PCASMSetTotalSubdomains(PC,int,IS *);
182 EXTERN int PCASMSetOverlap(PC,int);
183 /*E
184     PCASMType - Type of additive Schwarz method to use
185 
186 $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
187 $                 and computed values in ghost regions are added together. Classical
188 $                 standard additive Schwarz
189 $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
190 $                    region are discarded. Default
191 $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
192 $                       region are added back in
193 $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
194 $                not very good.
195 
196    Level: beginner
197 
198 .seealso: PCASMSetType()
199 E*/
200 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
201 EXTERN int PCASMSetType(PC,PCASMType);
202 EXTERN int PCASMCreateSubdomains2D(int,int,int,int,int,int,int *,IS **);
203 EXTERN int PCASMSetUseInPlace(PC);
204 EXTERN int PCASMGetLocalSubdomains(PC,int*,IS**);
205 
206 /*E
207     PCCompositeType - Determines how two or more preconditioner are composed
208 
209 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
210 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
211 $                                computed after the previous preconditioner application
212 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
213 $                         where first preconditioner is built from alpha I + S and second from
214 $                         alpha I + R
215 
216    Level: beginner
217 
218 .seealso: PCCompositeSetType()
219 E*/
220 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType;
221 EXTERN int PCCompositeSetType(PC,PCCompositeType);
222 EXTERN int PCCompositeAddPC(PC,PCType);
223 EXTERN int PCCompositeGetPC(PC pc,int n,PC *);
224 
225 EXTERN int PCRedundantSetScatter(PC,VecScatter,VecScatter);
226 EXTERN int PCRedundantGetOperators(PC,Mat*,Mat*);
227 EXTERN int PCRedundantGetPC(PC,PC*);
228 EXTERN int MatGetOrderingList(PetscFList *list);
229 #endif
230 
231 
232 
233 
234