Actual source code: petscpc.h
1: /*
2: Preconditioner module.
3: */
6: #include petscmat.h
9: EXTERN PetscErrorCode PCInitializePackage(const char[]);
11: /*
12: PCList contains the list of preconditioners currently registered
13: These are added with the PCRegisterDynamic() macro
14: */
17: /*S
18: PC - Abstract PETSc object that manages all preconditioners
20: Level: beginner
22: Concepts: preconditioners
24: .seealso: PCCreate(), PCSetType(), PCType (for list of available types)
25: S*/
26: typedef struct _p_PC* PC;
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()
33: Level: beginner
35: Notes: Click on the links below to see details on a particular solver
37: .seealso: PCSetType(), PC, PCCreate()
38: E*/
39: #define PCType const char*
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: #define PCGALERKIN "galerkin"
66: #define PCOPENMP "openmp"
67: #define PCSUPPORTGRAPH "supportgraph"
68: #define PCASA "asa"
70: /* Logging support */
73: /*E
74: PCSide - If the preconditioner is to be applied to the left, right
75: or symmetrically around the operator.
77: Level: beginner
79: .seealso:
80: E*/
81: typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
84: EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
85: EXTERN PetscErrorCode PCSetType(PC,PCType);
86: EXTERN PetscErrorCode PCSetUp(PC);
87: EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
88: EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
89: EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
90: EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
91: EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
92: EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
93: EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscTruth*);
94: EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
95: EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt);
96: EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscTruth*);
97: EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscTruth);
99: EXTERN PetscErrorCode PCRegisterDestroy(void);
100: EXTERN PetscErrorCode PCRegisterAll(const char[]);
103: EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
105: /*MC
106: PCRegisterDynamic - Adds a method to the preconditioner package.
108: Synopsis:
109: PetscErrorCode PCRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(PC))
111: Not collective
113: Input Parameters:
114: + name_solver - name of a new user-defined solver
115: . path - path (either absolute or relative) the library containing this solver
116: . name_create - name of routine to create method context
117: - routine_create - routine to create method context
119: Notes:
120: PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
122: If dynamic libraries are used, then the fourth input argument (routine_create)
123: is ignored.
125: Sample usage:
126: .vb
127: PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
128: "MySolverCreate",MySolverCreate);
129: .ve
131: Then, your solver can be chosen with the procedural interface via
132: $ PCSetType(pc,"my_solver")
133: or at runtime via the option
134: $ -pc_type my_solver
136: Level: advanced
138: Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, or ${any environmental variable}
139: occuring in pathname will be replaced with appropriate values.
140: If your function is not being put into a shared library then use PCRegister() instead
142: .keywords: PC, register
144: .seealso: PCRegisterAll(), PCRegisterDestroy()
145: M*/
146: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
147: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
148: #else
149: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
150: #endif
152: EXTERN PetscErrorCode PCDestroy(PC);
153: EXTERN PetscErrorCode PCSetFromOptions(PC);
154: EXTERN PetscErrorCode PCGetType(PC,PCType*);
156: EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
157: EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
158: EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
160: EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
161: EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
162: EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscTruth*,PetscTruth*);
164: EXTERN PetscErrorCode PCView(PC,PetscViewer);
166: EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
167: EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
168: EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);
170: EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
172: /*
173: These are used to provide extra scaling of preconditioned
174: operator for time-stepping schemes like in SUNDIALS
175: */
176: EXTERN PetscErrorCode PCDiagonalScale(PC,PetscTruth*);
177: EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
178: EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
179: EXTERN PetscErrorCode PCDiagonalScaleSet(PC,Vec);
181: /* ------------- options specific to particular preconditioners --------- */
183: EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
184: EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC);
185: EXTERN PetscErrorCode PCJacobiSetUseAbs(PC);
186: EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
187: EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
188: EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
190: EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
191: EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
193: #define USE_PRECONDITIONER_MATRIX 0
194: #define USE_TRUE_MATRIX 1
195: EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC);
196: EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
197: EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
199: EXTERN PetscErrorCode PCKSPSetUseTrue(PC);
201: EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(void*,Vec,Vec));
202: EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(void*,PCSide,Vec,Vec,Vec));
203: EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(void*,Vec,Vec));
204: EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(void*));
205: EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt));
206: EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(void*,PetscViewer));
207: EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(void*));
208: EXTERN PetscErrorCode PCShellGetContext(PC,void**);
209: EXTERN PetscErrorCode PCShellSetContext(PC,void*);
210: EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
211: EXTERN PetscErrorCode PCShellGetName(PC,char*[]);
213: EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);
214: EXTERN PetscErrorCode PCFactorSetShiftNonzero(PC,PetscReal);
215: EXTERN PetscErrorCode PCFactorSetShiftPd(PC,PetscTruth);
218: EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
219: EXTERN PetscErrorCode PCFactorSetPivoting(PC,PetscReal);
220: EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
222: EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType);
223: EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscTruth);
224: EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscTruth);
225: EXTERN PetscErrorCode PCFactorSetUseInPlace(PC);
226: EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
227: EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscTruth);
229: EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
230: EXTERN PetscErrorCode PCFactorSetUseDropTolerance(PC,PetscReal,PetscReal,PetscInt);
232: EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[]);
233: EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[]);
234: EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
235: /*E
236: PCASMType - Type of additive Schwarz method to use
238: $ PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
239: $ and computed values in ghost regions are added together. Classical
240: $ standard additive Schwarz
241: $ PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
242: $ region are discarded. Default
243: $ PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
244: $ region are added back in
245: $ PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
246: $ not very good.
248: Level: beginner
250: .seealso: PCASMSetType()
251: E*/
252: typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
255: EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
256: EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt *,IS **);
257: EXTERN PetscErrorCode PCASMSetUseInPlace(PC);
258: EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[]);
259: EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
261: /*E
262: PCCompositeType - Determines how two or more preconditioner are composed
264: $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
265: $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
266: $ computed after the previous preconditioner application
267: $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
268: $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
269: $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
270: $ where first preconditioner is built from alpha I + S and second from
271: $ alpha I + R
273: Level: beginner
275: .seealso: PCCompositeSetType()
276: E*/
277: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType;
280: EXTERN PetscErrorCode PCCompositeSetUseTrue(PC);
281: EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
282: EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
283: EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
284: EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
286: EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
287: EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
288: EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
289: EXTERN PetscErrorCode PCRedundantGetPC(PC,PC*);
291: EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
292: EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
293: EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
294: EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
295: EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
296: EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
297: EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
298: EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
300: EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
301: EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
302: EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
303: EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
305: EXTERN PetscErrorCode PCFieldSplitSetFields(PC,PetscInt,PetscInt*);
306: EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
307: EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
309: EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
310: EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
312: EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscReal*);
313: EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);
317: #endif /* __PETSCPC_H */