xref: /petsc/include/petscpc.h (revision 24c02a0f0f2ac65db2911f54fb0d656d2bff1739)
1d03aef70SBarry Smith /*
237f753daSBarry Smith       Preconditioner module.
3d03aef70SBarry Smith */
40a835dfdSSatish Balay #if !defined(__PETSCPC_H)
50a835dfdSSatish Balay #define __PETSCPC_H
60a835dfdSSatish Balay #include "petscmat.h"
7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
8d03aef70SBarry Smith 
9dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT  PCInitializePackage(const char[]);
101dbb0a54SBarry Smith 
11eec0b4cfSBarry Smith /*
12eec0b4cfSBarry Smith     PCList contains the list of preconditioners currently registered
13f1af5d2fSBarry Smith    These are added with the PCRegisterDynamic() macro
14eec0b4cfSBarry Smith */
15b0a32e0cSBarry Smith extern PetscFList PCList;
1682bf6240SBarry Smith 
173d957683SBarry Smith /*S
183d957683SBarry Smith      PC - Abstract PETSc object that manages all preconditioners
193d957683SBarry Smith 
203d957683SBarry Smith    Level: beginner
213d957683SBarry Smith 
223d957683SBarry Smith   Concepts: preconditioners
233d957683SBarry Smith 
241a480d89SAdministrator .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
253d957683SBarry Smith S*/
263d957683SBarry Smith typedef struct _p_PC* PC;
273d957683SBarry Smith 
283d957683SBarry Smith /*E
293d957683SBarry Smith     PCType - String with the name of a PETSc preconditioner method or the creation function
303d957683SBarry Smith        with an optional dynamic library name, for example
313d957683SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mypccreate()
323d957683SBarry Smith 
333d957683SBarry Smith    Level: beginner
343d957683SBarry Smith 
351a480d89SAdministrator    Notes: Click on the links below to see details on a particular solver
361a480d89SAdministrator 
371a480d89SAdministrator .seealso: PCSetType(), PC, PCCreate()
383d957683SBarry Smith E*/
39e5a9bf91SBarry Smith #define PCType const char*
4082bf6240SBarry Smith #define PCNONE            "none"
4182bf6240SBarry Smith #define PCJACOBI          "jacobi"
4282bf6240SBarry Smith #define PCSOR             "sor"
4382bf6240SBarry Smith #define PCLU              "lu"
4482bf6240SBarry Smith #define PCSHELL           "shell"
4582bf6240SBarry Smith #define PCBJACOBI         "bjacobi"
4682bf6240SBarry Smith #define PCMG              "mg"
4782bf6240SBarry Smith #define PCEISENSTAT       "eisenstat"
4882bf6240SBarry Smith #define PCILU             "ilu"
4982bf6240SBarry Smith #define PCICC             "icc"
5082bf6240SBarry Smith #define PCASM             "asm"
5194b7f48cSBarry Smith #define PCKSP             "ksp"
5282bf6240SBarry Smith #define PCCOMPOSITE       "composite"
53421c37bdSBarry Smith #define PCREDUNDANT       "redundant"
5427b520f0SBarry Smith #define PCSPAI            "spai"
55186905e3SBarry Smith #define PCNN              "nn"
564bbc92c1SBarry Smith #define PCCHOLESKY        "cholesky"
57a6b93e97SBarry Smith #define PCSAMG            "samg"
583a7fca6bSBarry Smith #define PCPBJACOBI        "pbjacobi"
597f5ff6fdSBarry Smith #define PCMAT             "mat"
60c4888f26SBarry Smith #define PCHYPRE           "hypre"
610971522cSBarry Smith #define PCFIELDSPLIT      "fieldsplit"
62be16f70fSSatish Balay #define PCTFS             "tfs"
635582bec1SHong Zhang #define PCML              "ml"
6436a49750SBarry Smith #define PCPROMETHEUS      "prometheus"
652a6744ebSBarry Smith #define PCGALERKIN        "galerkin"
66161c5ca9SBarry Smith #define PCOPENMP          "openmp"
67cf037197Ssdaitch #define PCSUPPORTGRAPH    "supportgraph"
68f4b8409dSBarry Smith #define PCASA             "asa"
69*24c02a0fSBarry Smith #define PCCP              "cp"
70123ea438SMatthew Knepley 
71123ea438SMatthew Knepley /* Logging support */
72dba47a55SKris Buschelman extern PetscCookie PETSCKSP_DLLEXPORT PC_COOKIE;
73123ea438SMatthew Knepley 
743d957683SBarry Smith /*E
753d957683SBarry Smith     PCSide - If the preconditioner is to be applied to the left, right
763d957683SBarry Smith      or symmetrically around the operator.
77d03aef70SBarry Smith 
783d957683SBarry Smith    Level: beginner
7921e95762SBarry Smith 
803d957683SBarry Smith .seealso:
813d957683SBarry Smith E*/
82aad2872bSLois Curfman McInnes typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
839dcbbd2bSBarry Smith extern const char *PCSides[];
8472b7852fSLois Curfman McInnes 
85dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCreate(MPI_Comm,PC*);
86f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetType(PC,PCType);
87dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetUp(PC);
88dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetUpOnBlocks(PC);
89dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApply(PC,Vec,Vec);
90dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricLeft(PC,Vec,Vec);
91dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricRight(PC,Vec,Vec);
92dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
93dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyTranspose(PC,Vec,Vec);
949cc050e5SLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyTransposeExists(PC,PetscTruth*);
95dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
96dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt);
97dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardsonExists(PC,PetscTruth*);
9891d066ecSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetInitialGuessNonzero(PC,PetscTruth);
9984cb2905SBarry Smith 
100dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterDestroy(void);
101dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterAll(const char[]);
1022bad1931SBarry Smith extern PetscTruth PCRegisterAllCalled;
10384cb2905SBarry Smith 
104dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
10530de9b25SBarry Smith 
10630de9b25SBarry Smith /*MC
10730de9b25SBarry Smith    PCRegisterDynamic - Adds a method to the preconditioner package.
10830de9b25SBarry Smith 
10930de9b25SBarry Smith    Synopsis:
110d360dc6fSBarry Smith    PetscErrorCode PCRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(PC))
11130de9b25SBarry Smith 
11230de9b25SBarry Smith    Not collective
11330de9b25SBarry Smith 
11430de9b25SBarry Smith    Input Parameters:
11530de9b25SBarry Smith +  name_solver - name of a new user-defined solver
11630de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
11730de9b25SBarry Smith .  name_create - name of routine to create method context
11830de9b25SBarry Smith -  routine_create - routine to create method context
11930de9b25SBarry Smith 
12030de9b25SBarry Smith    Notes:
12130de9b25SBarry Smith    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
12230de9b25SBarry Smith 
12330de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
12430de9b25SBarry Smith    is ignored.
12530de9b25SBarry Smith 
12630de9b25SBarry Smith    Sample usage:
12730de9b25SBarry Smith .vb
12830de9b25SBarry Smith    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
12930de9b25SBarry Smith               "MySolverCreate",MySolverCreate);
13030de9b25SBarry Smith .ve
13130de9b25SBarry Smith 
13230de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
13330de9b25SBarry Smith $     PCSetType(pc,"my_solver")
13430de9b25SBarry Smith    or at runtime via the option
13530de9b25SBarry Smith $     -pc_type my_solver
13630de9b25SBarry Smith 
13730de9b25SBarry Smith    Level: advanced
13830de9b25SBarry Smith 
139ab901514SBarry Smith    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
14030de9b25SBarry Smith            occuring in pathname will be replaced with appropriate values.
14130de9b25SBarry Smith          If your function is not being put into a shared library then use PCRegister() instead
14230de9b25SBarry Smith 
14330de9b25SBarry Smith .keywords: PC, register
14430de9b25SBarry Smith 
14530de9b25SBarry Smith .seealso: PCRegisterAll(), PCRegisterDestroy()
14630de9b25SBarry Smith M*/
147aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
148f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
1496df38c32SLois Curfman McInnes #else
150f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
1516df38c32SLois Curfman McInnes #endif
1526df38c32SLois Curfman McInnes 
153dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDestroy(PC);
154dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetFromOptions(PC);
155dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetType(PC,PCType*);
15614c91fddSBarry Smith 
157a4fd02acSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatrix(PC,Mat*);
158dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
159dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
1605b116368SBarry Smith 
161dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetOperators(PC,Mat,Mat,MatStructure);
162dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperators(PC,Mat*,Mat*,MatStructure*);
163906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperatorsSet(PC,PetscTruth*,PetscTruth*);
1644b0e389bSBarry Smith 
165dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCView(PC,PetscViewer);
1667bc3d0afSSatish Balay 
167dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetOptionsPrefix(PC,const char[]);
168dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCAppendOptionsPrefix(PC,const char[]);
1692dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOptionsPrefix(PC,const char*[]);
1708ed539a5SBarry Smith 
171dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCComputeExplicitOperator(PC,Mat*);
17271601f6fSBarry Smith 
173d6913704SBarry Smith /*
174d6913704SBarry Smith       These are used to provide extra scaling of preconditioned
1750f3b3ca1SHong Zhang    operator for time-stepping schemes like in SUNDIALS
176d6913704SBarry Smith */
177dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScale(PC,PetscTruth*);
178dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleLeft(PC,Vec,Vec);
179dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleRight(PC,Vec,Vec);
180dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleSet(PC,Vec);
181d6913704SBarry Smith 
18284cb2905SBarry Smith /* ------------- options specific to particular preconditioners --------- */
183329f5518SBarry Smith 
184dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowMax(PC);
185d87ca0b3SLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowSum(PC);
186cd47f5d9SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseAbs(PC);
187dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric(PC,MatSORType);
188dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega(PC,PetscReal);
189dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations(PC,PetscInt,PetscInt);
190d03aef70SBarry Smith 
191dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatSetOmega(PC,PetscReal);
192dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatNoDiagonalScaling(PC);
193421c37bdSBarry Smith 
19415aa81f8SBarry Smith #define USE_PRECONDITIONER_MATRIX 0
19515aa81f8SBarry Smith #define USE_TRUE_MATRIX           1
196dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetUseTrueLocal(PC);
197dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
198dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
1991eb62cbbSBarry Smith 
200dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPSetUseTrue(PC);
201981c4779SBarry Smith 
202be29d3c6SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApply(PC,PetscErrorCode (*)(void*,Vec,Vec));
2032bb17772SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyBA(PC,PetscErrorCode (*)(void*,PCSide,Vec,Vec,Vec));
204dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyTranspose(PC,PetscErrorCode (*)(void*,Vec,Vec));
205dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetSetUp(PC,PetscErrorCode (*)(void*));
206be29d3c6SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyRichardson(PC,PetscErrorCode (*)(void*,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt));
207dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetView(PC,PetscErrorCode (*)(void*,PetscViewer));
20818be62a5SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetDestroy(PC,PetscErrorCode (*)(void*));
209be29d3c6SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetContext(PC,void**);
210be29d3c6SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetContext(PC,void*);
211dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetName(PC,const char[]);
212dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetName(PC,char*[]);
213aabeff55SBarry Smith 
214dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot(PC,PetscReal);
215dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero(PC,PetscReal);
216dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd(PC,PetscTruth);
217ee45ca4aSHong Zhang 
2182401956bSBarry Smith 
21955ba2a51SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC,PetscReal);
2202401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC,PetscReal);
2212401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
222421c37bdSBarry Smith 
223e5a9bf91SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC,MatOrderingType);
2242401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC,PetscTruth);
2252401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC,PetscTruth);
2262401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC);
22712a54f56SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC);
2282401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC,PetscTruth);
229f5a88c2aSLois Curfman McInnes 
2302401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC,PetscInt);
2312401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC,PetscReal,PetscReal,PetscInt);
232b35a507dSBarry Smith 
233dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetLocalSubdomains(PC,PetscInt,IS[]);
234dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetTotalSubdomains(PC,PetscInt,IS[]);
235dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetOverlap(PC,PetscInt);
2363d957683SBarry Smith /*E
2373d957683SBarry Smith     PCASMType - Type of additive Schwarz method to use
2383d957683SBarry Smith 
2393d957683SBarry Smith $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
2403d957683SBarry Smith $                 and computed values in ghost regions are added together. Classical
2413d957683SBarry Smith $                 standard additive Schwarz
2423d957683SBarry Smith $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
2433d957683SBarry Smith $                    region are discarded. Default
2443d957683SBarry Smith $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
2453d957683SBarry Smith $                       region are added back in
2463d957683SBarry Smith $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
2473d957683SBarry Smith $                not very good.
2483d957683SBarry Smith 
2493d957683SBarry Smith    Level: beginner
2503d957683SBarry Smith 
2513d957683SBarry Smith .seealso: PCASMSetType()
2523d957683SBarry Smith E*/
253d252947aSBarry Smith typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
2549dcbbd2bSBarry Smith extern const char *PCASMTypes[];
2559dcbbd2bSBarry Smith 
256dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetType(PC,PCASMType);
257dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt *,IS **);
258dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetUseInPlace(PC);
259dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubdomains(PC,PetscInt*,IS*[]);
260dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
261981c4779SBarry Smith 
2623d957683SBarry Smith /*E
2633d957683SBarry Smith     PCCompositeType - Determines how two or more preconditioner are composed
2643d957683SBarry Smith 
2653d957683SBarry Smith $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
2663d957683SBarry Smith $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
2673d957683SBarry Smith $                                computed after the previous preconditioner application
268ccb205f8SBarry Smith $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
269ccb205f8SBarry Smith $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
2703d957683SBarry Smith $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
2713d957683SBarry Smith $                         where first preconditioner is built from alpha I + S and second from
2723d957683SBarry Smith $                         alpha I + R
2733d957683SBarry Smith 
2743d957683SBarry Smith    Level: beginner
2753d957683SBarry Smith 
2763d957683SBarry Smith .seealso: PCCompositeSetType()
2773d957683SBarry Smith E*/
27851f519a2SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL} PCCompositeType;
2799dcbbd2bSBarry Smith extern const char *PCCompositeTypes[];
2809dcbbd2bSBarry Smith 
281dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC);
282dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC,PCCompositeType);
283dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC,PCType);
2849a05e725SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC,PetscInt,PC *);
285dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC,PetscScalar);
286981c4779SBarry Smith 
28709a6bc64SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber(PC,PetscInt);
288dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter(PC,VecScatter,VecScatter);
289dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators(PC,Mat*,Mat*);
290dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC(PC,PC*);
291da3a660dSBarry Smith 
292dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetEpsilon(PC,double);
293dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetNBSteps(PC,PetscInt);
294dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMax(PC,PetscInt);
295dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMaxNew(PC,PetscInt);
296dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetBlockSize(PC,PetscInt);
297dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetCacheSize(PC,PetscInt);
298dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetVerbose(PC,PetscInt);
299dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetSp(PC,PetscInt);
3003304466cSBarry Smith 
301dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC,const char[]);
302dd896f7dSSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC,const char*[]);
303dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
304dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
3053304466cSBarry Smith 
306dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC,PetscInt,PetscInt*);
307dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC,PCCompositeType);
30851f519a2SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC,PetscInt);
3093d30b1ffSBarry Smith 
3102a6744ebSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetRestriction(PC,Mat);
3112a6744ebSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetInterpolation(PC,Mat);
3122a6744ebSBarry Smith 
3131a3fd43dSSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetCoordinates(PC,PetscInt,PetscReal*);
3142102ba4dSSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSASetVectors(PC,PetscInt,PetscReal *);
3152102ba4dSSatish Balay 
31667fac13cSBarry Smith 
317e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
318123ea438SMatthew Knepley #endif /* __PETSCPC_H */
319