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"
68: /* Logging support */
71: /*E
72: PCSide - If the preconditioner is to be applied to the left, right
73: or symmetrically around the operator.
75: Level: beginner
77: .seealso:
78: E*/
79: typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
82: EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
83: EXTERN PetscErrorCode PCSetType(PC,PCType);
84: EXTERN PetscErrorCode PCSetUp(PC);
85: EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
86: EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
87: EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
88: EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
89: EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
90: EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
91: EXTERN PetscErrorCode PCHasApplyTranspose(PC,PetscTruth*);
92: EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
93: EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt);
94: EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscTruth*);
96: EXTERN PetscErrorCode PCRegisterDestroy(void);
97: EXTERN PetscErrorCode PCRegisterAll(const char[]);
100: EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
102: /*MC
103: PCRegisterDynamic - Adds a method to the preconditioner package.
105: Synopsis:
106: PetscErrorCode PCRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(PC))
108: Not collective
110: Input Parameters:
111: + name_solver - name of a new user-defined solver
112: . path - path (either absolute or relative) the library containing this solver
113: . name_create - name of routine to create method context
114: - routine_create - routine to create method context
116: Notes:
117: PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
119: If dynamic libraries are used, then the fourth input argument (routine_create)
120: is ignored.
122: Sample usage:
123: .vb
124: PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
125: "MySolverCreate",MySolverCreate);
126: .ve
128: Then, your solver can be chosen with the procedural interface via
129: $ PCSetType(pc,"my_solver")
130: or at runtime via the option
131: $ -pc_type my_solver
133: Level: advanced
135: Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, or ${any environmental variable}
136: occuring in pathname will be replaced with appropriate values.
137: If your function is not being put into a shared library then use PCRegister() instead
139: .keywords: PC, register
141: .seealso: PCRegisterAll(), PCRegisterDestroy()
142: M*/
143: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
144: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
145: #else
146: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
147: #endif
149: EXTERN PetscErrorCode PCDestroy(PC);
150: EXTERN PetscErrorCode PCSetFromOptions(PC);
151: EXTERN PetscErrorCode PCGetType(PC,PCType*);
153: EXTERN PetscErrorCode PCGetFactoredMatrix(PC,Mat*);
154: EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
155: EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
157: EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
158: EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
159: EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscTruth*,PetscTruth*);
161: EXTERN PetscErrorCode PCView(PC,PetscViewer);
163: EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
164: EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
165: EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);
167: EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
169: /*
170: These are used to provide extra scaling of preconditioned
171: operator for time-stepping schemes like in SUNDIALS
172: */
173: EXTERN PetscErrorCode PCDiagonalScale(PC,PetscTruth*);
174: EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
175: EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
176: EXTERN PetscErrorCode PCDiagonalScaleSet(PC,Vec);
178: /* ------------- options specific to particular preconditioners --------- */
180: EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
181: EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC);
182: EXTERN PetscErrorCode PCJacobiSetUseAbs(PC);
183: EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
184: EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
185: EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
187: EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
188: EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
190: #define USE_PRECONDITIONER_MATRIX 0
191: #define USE_TRUE_MATRIX 1
192: EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC);
193: EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
194: EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
196: EXTERN PetscErrorCode PCKSPSetUseTrue(PC);
198: EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(void*,Vec,Vec));
199: EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(void*,PCSide,Vec,Vec,Vec));
200: EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(void*,Vec,Vec));
201: EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(void*));
202: EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt));
203: EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(void*,PetscViewer));
204: EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(void*));
205: EXTERN PetscErrorCode PCShellGetContext(PC,void**);
206: EXTERN PetscErrorCode PCShellSetContext(PC,void*);
207: EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
208: EXTERN PetscErrorCode PCShellGetName(PC,char*[]);
210: EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);
211: EXTERN PetscErrorCode PCFactorSetShiftNonzero(PC,PetscReal);
212: EXTERN PetscErrorCode PCFactorSetShiftPd(PC,PetscTruth);
215: EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
216: EXTERN PetscErrorCode PCFactorSetPivoting(PC,PetscReal);
217: EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
219: EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType);
220: EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscTruth);
221: EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscTruth);
222: EXTERN PetscErrorCode PCFactorSetUseInPlace(PC);
223: EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
224: EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscTruth);
226: EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
227: EXTERN PetscErrorCode PCFactorSetUseDropTolerance(PC,PetscReal,PetscReal,PetscInt);
229: EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[]);
230: EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[]);
231: EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
232: /*E
233: PCASMType - Type of additive Schwarz method to use
235: $ PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
236: $ and computed values in ghost regions are added together. Classical
237: $ standard additive Schwarz
238: $ PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
239: $ region are discarded. Default
240: $ PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
241: $ region are added back in
242: $ PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
243: $ not very good.
245: Level: beginner
247: .seealso: PCASMSetType()
248: E*/
249: typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
252: EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
253: EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt *,IS **);
254: EXTERN PetscErrorCode PCASMSetUseInPlace(PC);
255: EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[]);
256: EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
258: /*E
259: PCCompositeType - Determines how two or more preconditioner are composed
261: $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
262: $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
263: $ computed after the previous preconditioner application
264: $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
265: $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
266: $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
267: $ where first preconditioner is built from alpha I + S and second from
268: $ alpha I + R
270: Level: beginner
272: .seealso: PCCompositeSetType()
273: E*/
274: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType;
277: EXTERN PetscErrorCode PCCompositeSetUseTrue(PC);
278: EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
279: EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
280: EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
281: EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
283: EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
284: EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
285: EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
286: EXTERN PetscErrorCode PCRedundantGetPC(PC,PC*);
288: EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
289: EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
290: EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
291: EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
292: EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
293: EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
294: EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
295: EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
297: EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
298: EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
299: EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
300: EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
302: EXTERN PetscErrorCode PCFieldSplitSetFields(PC,PetscInt,PetscInt*);
303: EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
304: EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
306: EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
307: EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
309: EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscReal*);
310: EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);
314: #endif /* __PETSCPC_H */