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 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 PCEXOTIC "exotic"
67: #define PCOPENMP "openmp"
68: #define PCSUPPORTGRAPH "supportgraph"
69: #define PCASA "asa"
70: #define PCCP "cp"
71: #define PCBFBT "bfbt"
72: #define PCPYTHON "python"
74: /* Logging support */
77: /*E
78: PCSide - If the preconditioner is to be applied to the left, right
79: or symmetrically around the operator.
81: Level: beginner
83: .seealso:
84: E*/
85: typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
88: EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
89: EXTERN PetscErrorCode PCSetType(PC,const PCType);
90: EXTERN PetscErrorCode PCSetUp(PC);
91: EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
92: EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
93: EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
94: EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
95: EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
96: EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
97: EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscTruth*);
98: EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
100: /*E
101: PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
103: Level: advanced
105: Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
107: .seealso: PCApplyRichardson()
108: E*/
109: typedef enum {
110: PCRICHARDSON_CONVERGED_RTOL = 2,
111: PCRICHARDSON_CONVERGED_ATOL = 3,
112: PCRICHARDSON_CONVERGED_ITS = 4,
113: PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason;
115: EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscInt*,PCRichardsonConvergedReason*);
116: EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscTruth*);
117: EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscTruth);
119: EXTERN PetscErrorCode PCRegisterDestroy(void);
120: EXTERN PetscErrorCode PCRegisterAll(const char[]);
123: EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
125: /*MC
126: PCRegisterDynamic - Adds a method to the preconditioner package.
128: Synopsis:
129: PetscErrorCode PCRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(PC))
131: Not collective
133: Input Parameters:
134: + name_solver - name of a new user-defined solver
135: . path - path (either absolute or relative) the library containing this solver
136: . name_create - name of routine to create method context
137: - routine_create - routine to create method context
139: Notes:
140: PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
142: If dynamic libraries are used, then the fourth input argument (routine_create)
143: is ignored.
145: Sample usage:
146: .vb
147: PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
148: "MySolverCreate",MySolverCreate);
149: .ve
151: Then, your solver can be chosen with the procedural interface via
152: $ PCSetType(pc,"my_solver")
153: or at runtime via the option
154: $ -pc_type my_solver
156: Level: advanced
158: Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, or ${any environmental variable}
159: occuring in pathname will be replaced with appropriate values.
160: If your function is not being put into a shared library then use PCRegister() instead
162: .keywords: PC, register
164: .seealso: PCRegisterAll(), PCRegisterDestroy()
165: M*/
166: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
167: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
168: #else
169: #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
170: #endif
172: EXTERN PetscErrorCode PCDestroy(PC);
173: EXTERN PetscErrorCode PCSetFromOptions(PC);
174: EXTERN PetscErrorCode PCGetType(PC,const PCType*);
176: EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
177: EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
178: EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
180: EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
181: EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
182: EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscTruth*,PetscTruth*);
184: EXTERN PetscErrorCode PCView(PC,PetscViewer);
186: EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
187: EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
188: EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);
190: EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
192: /*
193: These are used to provide extra scaling of preconditioned
194: operator for time-stepping schemes like in SUNDIALS
195: */
196: EXTERN PetscErrorCode PCDiagonalScale(PC,PetscTruth*);
197: EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
198: EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
199: EXTERN PetscErrorCode PCDiagonalScaleSet(PC,Vec);
201: /* ------------- options specific to particular preconditioners --------- */
203: EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
204: EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC);
205: EXTERN PetscErrorCode PCJacobiSetUseAbs(PC);
206: EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
207: EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
208: EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
210: EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
211: EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
213: #define USE_PRECONDITIONER_MATRIX 0
214: #define USE_TRUE_MATRIX 1
215: EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC);
216: EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
217: EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
219: EXTERN PetscErrorCode PCKSPSetUseTrue(PC);
221: EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(void*,Vec,Vec));
222: EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(void*,PCSide,Vec,Vec,Vec));
223: EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(void*,Vec,Vec));
224: EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(void*));
225: EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscInt*,PCRichardsonConvergedReason*));
226: EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(void*,PetscViewer));
227: EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(void*));
228: EXTERN PetscErrorCode PCShellGetContext(PC,void**);
229: EXTERN PetscErrorCode PCShellSetContext(PC,void*);
230: EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
231: EXTERN PetscErrorCode PCShellGetName(PC,char*[]);
233: EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);
234: EXTERN PetscErrorCode PCFactorSetShiftNonzero(PC,PetscReal);
235: EXTERN PetscErrorCode PCFactorSetShiftPd(PC,PetscTruth);
236: EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
237: EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
239: EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
240: EXTERN PetscErrorCode PCFactorSetPivoting(PC,PetscReal);
241: EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
243: EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,const MatOrderingType);
244: EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscTruth);
245: EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscTruth);
246: EXTERN PetscErrorCode PCFactorSetUseInPlace(PC);
247: EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
248: EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscTruth);
249: EXTERN PetscErrorCode PCFactorSetShiftInBlocks(PC,PetscReal);
251: EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
252: EXTERN PetscErrorCode PCFactorSetUseDropTolerance(PC,PetscReal,PetscReal,PetscInt);
254: EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[]);
255: EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[]);
256: EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
257: /*E
258: PCASMType - Type of additive Schwarz method to use
260: $ PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
261: $ and computed values in ghost regions are added together. Classical
262: $ standard additive Schwarz
263: $ PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
264: $ region are discarded. Default
265: $ PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
266: $ region are added back in
267: $ PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
268: $ not very good.
270: Level: beginner
272: .seealso: PCASMSetType()
273: E*/
274: typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
277: EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
278: EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
279: EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[]);
280: EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt *,IS **);
281: EXTERN PetscErrorCode PCASMSetUseInPlace(PC);
282: EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[]);
283: EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
285: /*E
286: PCCompositeType - Determines how two or more preconditioner are composed
288: $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
289: $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
290: $ computed after the previous preconditioner application
291: $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
292: $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
293: $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
294: $ where first preconditioner is built from alpha I + S and second from
295: $ alpha I + R
297: Level: beginner
299: .seealso: PCCompositeSetType()
300: E*/
301: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
304: EXTERN PetscErrorCode PCCompositeSetUseTrue(PC);
305: EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
306: EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
307: EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
308: EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
310: EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
311: EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
312: EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
313: EXTERN PetscErrorCode PCRedundantGetPC(PC,PC*);
315: EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
316: EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
317: EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
318: EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
319: EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
320: EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
321: EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
322: EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
324: EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
325: EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
326: EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
327: EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
329: EXTERN PetscErrorCode PCFieldSplitSetFields(PC,PetscInt,PetscInt*);
330: EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
331: EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
332: EXTERN PetscErrorCode PCFieldSplitSetIS(PC,IS);
333: EXTERN PetscErrorCode PCFieldSplitSchurPrecondition(PC,PetscTruth);
334: EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
336: EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
337: EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
339: EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscReal*);
340: EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);
342: EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);
345: #endif /* __PETSCPC_H */