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: */
16: #define PCType const char*
18: /*S
19: PC - Abstract PETSc object that manages all preconditioners
21: Level: beginner
23: Concepts: preconditioners
25: .seealso: PCCreate(), PCSetType(), PCType (for list of available types)
26: S*/
27: typedef struct _p_PC* PC;
29: /*E
30: PCType - String with the name of a PETSc preconditioner method or the creation function
31: with an optional dynamic library name, for example
32: http://www.mcs.anl.gov/petsc/lib.a:mypccreate()
34: Level: beginner
36: Notes: Click on the links below to see details on a particular solver
38: .seealso: PCSetType(), PC, PCCreate()
39: E*/
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"
67: /* Logging support */
70: /*E
71: PCSide - If the preconditioner is to be applied to the left, right
72: or symmetrically around the operator.
74: Level: beginner
76: .seealso:
77: E*/
78: typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
81: EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
82: EXTERN PetscErrorCode PCSetType(PC,PCType);
83: EXTERN PetscErrorCode PCSetUp(PC);
84: EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
85: EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
86: EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
87: EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
88: EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
89: EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
90: EXTERN PetscErrorCode PCHasApplyTranspose(PC,PetscTruth*);
91: EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
92: EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt);
93: EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscTruth*);
95: EXTERN PetscErrorCode PCRegisterDestroy(void);
96: EXTERN PetscErrorCode PCRegisterAll(const char[]);
99: EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
101: /*MC
102: PCRegisterDynamic - Adds a method to the preconditioner package.
104: Synopsis:
105: PetscErrorCode PCRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(PC))
107: Not collective
109: Input Parameters:
110: + name_solver - name of a new user-defined solver
111: . path - path (either absolute or relative) the library containing this solver
112: . name_create - name of routine to create method context
113: - routine_create - routine to create method context
115: Notes:
116: PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
118: If dynamic libraries are used, then the fourth input argument (routine_create)
119: is ignored.
121: Sample usage:
122: .vb
123: PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
124: "MySolverCreate",MySolverCreate);
125: .ve
127: Then, your solver can be chosen with the procedural interface via
128: $ PCSetType(pc,"my_solver")
129: or at runtime via the option
130: $ -pc_type my_solver
132: Level: advanced
134: Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, or ${any environmental variable}
135: occuring in pathname will be replaced with appropriate values.
136: If your function is not being put into a shared library then use PCRegister() instead
138: .keywords: PC, register
140: .seealso: PCRegisterAll(), PCRegisterDestroy()
141: M*/
142: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
143: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
144: #else
145: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
146: #endif
148: EXTERN PetscErrorCode PCDestroy(PC);
149: EXTERN PetscErrorCode PCSetFromOptions(PC);
150: EXTERN PetscErrorCode PCGetType(PC,PCType*);
152: EXTERN PetscErrorCode PCGetFactoredMatrix(PC,Mat*);
153: EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
154: EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
156: EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
157: EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
158: EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscTruth*,PetscTruth*);
160: EXTERN PetscErrorCode PCView(PC,PetscViewer);
162: EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
163: EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
164: EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);
166: EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
168: /*
169: These are used to provide extra scaling of preconditioned
170: operator for time-stepping schemes like in SUNDIALS
171: */
172: EXTERN PetscErrorCode PCDiagonalScale(PC,PetscTruth*);
173: EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
174: EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
175: EXTERN PetscErrorCode PCDiagonalScaleSet(PC,Vec);
177: /* ------------- options specific to particular preconditioners --------- */
179: EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
180: EXTERN PetscErrorCode PCJacobiSetUseAbs(PC);
181: EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
182: EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
183: EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
185: EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
186: EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
188: #define USE_PRECONDITIONER_MATRIX 0
189: #define USE_TRUE_MATRIX 1
190: EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC);
191: EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
192: EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
194: EXTERN PetscErrorCode PCKSPSetUseTrue(PC);
196: EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(void*,Vec,Vec));
197: EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(void*,PCSide,Vec,Vec,Vec));
198: EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(void*,Vec,Vec));
199: EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(void*));
200: EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt));
201: EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(void*,PetscViewer));
202: EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(void*));
203: EXTERN PetscErrorCode PCShellGetContext(PC,void**);
204: EXTERN PetscErrorCode PCShellSetContext(PC,void*);
205: EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
206: EXTERN PetscErrorCode PCShellGetName(PC,char*[]);
208: EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);
209: EXTERN PetscErrorCode PCFactorSetShiftNonzero(PC,PetscReal);
210: EXTERN PetscErrorCode PCFactorSetShiftPd(PC,PetscTruth);
213: EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
214: EXTERN PetscErrorCode PCFactorSetPivoting(PC,PetscReal);
215: EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
217: EXTERN PetscErrorCode PCFactorSetMatOrdering(PC,MatOrderingType);
218: EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscTruth);
219: EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscTruth);
220: EXTERN PetscErrorCode PCFactorSetUseInPlace(PC);
221: EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
222: EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscTruth);
224: EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
225: EXTERN PetscErrorCode PCFactorSetUseDropTolerance(PC,PetscReal,PetscReal,PetscInt);
227: EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[]);
228: EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[]);
229: EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
230: /*E
231: PCASMType - Type of additive Schwarz method to use
233: $ PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
234: $ and computed values in ghost regions are added together. Classical
235: $ standard additive Schwarz
236: $ PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
237: $ region are discarded. Default
238: $ PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
239: $ region are added back in
240: $ PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
241: $ not very good.
243: Level: beginner
245: .seealso: PCASMSetType()
246: E*/
247: typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
250: EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
251: EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt *,IS **);
252: EXTERN PetscErrorCode PCASMSetUseInPlace(PC);
253: EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[]);
254: EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
256: /*E
257: PCCompositeType - Determines how two or more preconditioner are composed
259: $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
260: $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
261: $ computed after the previous preconditioner application
262: $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
263: $ where first preconditioner is built from alpha I + S and second from
264: $ alpha I + R
266: Level: beginner
268: .seealso: PCCompositeSetType()
269: E*/
270: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType;
273: EXTERN PetscErrorCode PCCompositeSetUseTrue(PC);
274: EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
275: EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
276: EXTERN PetscErrorCode PCCompositeGetPC(PC pc,PetscInt n,PC *);
277: EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
279: EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
280: EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
281: EXTERN PetscErrorCode PCRedundantGetPC(PC,PC*);
283: EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
284: EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
285: EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
286: EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
287: EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
288: EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
289: EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
290: EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
292: EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
293: EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
294: EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
295: EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
297: EXTERN PetscErrorCode PCFieldSplitSetFields(PC,PetscInt,PetscInt*);
298: EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
300: EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
301: EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
303: EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscReal*);
304: EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);
308: #endif /* __PETSCPC_H */