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