xref: /petsc/include/petscpc.h (revision 0700a8246d308f50502909ba325e6169d3ee27eb)
1d03aef70SBarry Smith /*
237f753daSBarry Smith       Preconditioner module.
3d03aef70SBarry Smith */
40a835dfdSSatish Balay #if !defined(__PETSCPC_H)
50a835dfdSSatish Balay #define __PETSCPC_H
60a835dfdSSatish Balay #include "petscmat.h"
76c699258SBarry Smith #include "petscda.h"
8e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
9d03aef70SBarry Smith 
10dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT  PCInitializePackage(const char[]);
111dbb0a54SBarry Smith 
12eec0b4cfSBarry Smith /*
13eec0b4cfSBarry Smith     PCList contains the list of preconditioners currently registered
14f1af5d2fSBarry Smith    These are added with the PCRegisterDynamic() macro
15eec0b4cfSBarry Smith */
16b0a32e0cSBarry Smith extern PetscFList PCList;
1782bf6240SBarry Smith 
183d957683SBarry Smith /*S
193d957683SBarry Smith      PC - Abstract PETSc object that manages all preconditioners
203d957683SBarry Smith 
213d957683SBarry Smith    Level: beginner
223d957683SBarry Smith 
233d957683SBarry Smith   Concepts: preconditioners
243d957683SBarry Smith 
251a480d89SAdministrator .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
263d957683SBarry Smith S*/
273d957683SBarry Smith typedef struct _p_PC* PC;
283d957683SBarry Smith 
293d957683SBarry Smith /*E
303d957683SBarry Smith     PCType - String with the name of a PETSc preconditioner method or the creation function
313d957683SBarry Smith        with an optional dynamic library name, for example
323d957683SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mypccreate()
333d957683SBarry Smith 
343d957683SBarry Smith    Level: beginner
353d957683SBarry Smith 
361a480d89SAdministrator    Notes: Click on the links below to see details on a particular solver
371a480d89SAdministrator 
381a480d89SAdministrator .seealso: PCSetType(), PC, PCCreate()
393d957683SBarry Smith E*/
40a313700dSBarry Smith #define PCType char*
4182bf6240SBarry Smith #define PCNONE            "none"
4282bf6240SBarry Smith #define PCJACOBI          "jacobi"
4382bf6240SBarry Smith #define PCSOR             "sor"
4482bf6240SBarry Smith #define PCLU              "lu"
4582bf6240SBarry Smith #define PCSHELL           "shell"
4682bf6240SBarry Smith #define PCBJACOBI         "bjacobi"
4782bf6240SBarry Smith #define PCMG              "mg"
4882bf6240SBarry Smith #define PCEISENSTAT       "eisenstat"
4982bf6240SBarry Smith #define PCILU             "ilu"
5082bf6240SBarry Smith #define PCICC             "icc"
5182bf6240SBarry Smith #define PCASM             "asm"
5294b7f48cSBarry Smith #define PCKSP             "ksp"
5382bf6240SBarry Smith #define PCCOMPOSITE       "composite"
54421c37bdSBarry Smith #define PCREDUNDANT       "redundant"
5527b520f0SBarry Smith #define PCSPAI            "spai"
56186905e3SBarry Smith #define PCNN              "nn"
574bbc92c1SBarry Smith #define PCCHOLESKY        "cholesky"
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"
667233f9f0SBarry Smith #define PCEXOTIC          "exotic"
67161c5ca9SBarry Smith #define PCOPENMP          "openmp"
68cf037197Ssdaitch #define PCSUPPORTGRAPH    "supportgraph"
69f4b8409dSBarry Smith #define PCASA             "asa"
7024c02a0fSBarry Smith #define PCCP              "cp"
71628b8437SMatthew Knepley #define PCBFBT            "bfbt"
72519d70e2SJed Brown #define PCLSC             "lsc"
731d6018f0SLisandro Dalcin #define PCPYTHON          "python"
74f91d8e95SBarry Smith #define PCPFMG            "pfmg"
75d851a50bSGlenn Hammond #define PCSYSPFMG         "syspfmg"
76df826632SBarry Smith #define PCREDISTRIBUTE    "redistribute"
77123ea438SMatthew Knepley 
78123ea438SMatthew Knepley /* Logging support */
79*0700a824SBarry Smith extern PetscClassId PETSCKSP_DLLEXPORT PC_CLASSID;
80123ea438SMatthew Knepley 
813d957683SBarry Smith /*E
823d957683SBarry Smith     PCSide - If the preconditioner is to be applied to the left, right
833d957683SBarry Smith      or symmetrically around the operator.
84d03aef70SBarry Smith 
853d957683SBarry Smith    Level: beginner
8621e95762SBarry Smith 
873d957683SBarry Smith .seealso:
883d957683SBarry Smith E*/
89aad2872bSLois Curfman McInnes typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
909dcbbd2bSBarry Smith extern const char *PCSides[];
9172b7852fSLois Curfman McInnes 
92dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCreate(MPI_Comm,PC*);
93a313700dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetType(PC,const PCType);
94dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetUp(PC);
95dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetUpOnBlocks(PC);
96dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApply(PC,Vec,Vec);
97dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricLeft(PC,Vec,Vec);
98dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplySymmetricRight(PC,Vec,Vec);
99dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
100dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyTranspose(PC,Vec,Vec);
1019cc050e5SLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyTransposeExists(PC,PetscTruth*);
102dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
1034d0a8057SBarry Smith 
1044d0a8057SBarry Smith /*E
1054d0a8057SBarry Smith     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
1064d0a8057SBarry Smith 
1074d0a8057SBarry Smith    Level: advanced
1084d0a8057SBarry Smith 
1094d0a8057SBarry Smith    Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
1104d0a8057SBarry Smith 
1114d0a8057SBarry Smith .seealso: PCApplyRichardson()
1124d0a8057SBarry Smith E*/
1134d0a8057SBarry Smith typedef enum {
1144d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_RTOL               =  2,
1154d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_ATOL               =  3,
1164d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_ITS                =  4,
1174d0a8057SBarry Smith               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
1184d0a8057SBarry Smith 
1197319c654SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscTruth,PetscInt*,PCRichardsonConvergedReason*);
120dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCApplyRichardsonExists(PC,PetscTruth*);
12191d066ecSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetInitialGuessNonzero(PC,PetscTruth);
12284cb2905SBarry Smith 
123dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterDestroy(void);
124dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegisterAll(const char[]);
1252bad1931SBarry Smith extern PetscTruth PCRegisterAllCalled;
12684cb2905SBarry Smith 
127dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
12830de9b25SBarry Smith 
12930de9b25SBarry Smith /*MC
13030de9b25SBarry Smith    PCRegisterDynamic - Adds a method to the preconditioner package.
13130de9b25SBarry Smith 
13230de9b25SBarry Smith    Synopsis:
133d360dc6fSBarry Smith    PetscErrorCode PCRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(PC))
13430de9b25SBarry Smith 
13530de9b25SBarry Smith    Not collective
13630de9b25SBarry Smith 
13730de9b25SBarry Smith    Input Parameters:
13830de9b25SBarry Smith +  name_solver - name of a new user-defined solver
13930de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
14030de9b25SBarry Smith .  name_create - name of routine to create method context
14130de9b25SBarry Smith -  routine_create - routine to create method context
14230de9b25SBarry Smith 
14330de9b25SBarry Smith    Notes:
14430de9b25SBarry Smith    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
14530de9b25SBarry Smith 
14630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
14730de9b25SBarry Smith    is ignored.
14830de9b25SBarry Smith 
14930de9b25SBarry Smith    Sample usage:
15030de9b25SBarry Smith .vb
15130de9b25SBarry Smith    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
15230de9b25SBarry Smith               "MySolverCreate",MySolverCreate);
15330de9b25SBarry Smith .ve
15430de9b25SBarry Smith 
15530de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
15630de9b25SBarry Smith $     PCSetType(pc,"my_solver")
15730de9b25SBarry Smith    or at runtime via the option
15830de9b25SBarry Smith $     -pc_type my_solver
15930de9b25SBarry Smith 
16030de9b25SBarry Smith    Level: advanced
16130de9b25SBarry Smith 
162ab901514SBarry Smith    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
16330de9b25SBarry Smith            occuring in pathname will be replaced with appropriate values.
16430de9b25SBarry Smith          If your function is not being put into a shared library then use PCRegister() instead
16530de9b25SBarry Smith 
16630de9b25SBarry Smith .keywords: PC, register
16730de9b25SBarry Smith 
16830de9b25SBarry Smith .seealso: PCRegisterAll(), PCRegisterDestroy()
16930de9b25SBarry Smith M*/
170aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
171f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
1726df38c32SLois Curfman McInnes #else
173f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
1746df38c32SLois Curfman McInnes #endif
1756df38c32SLois Curfman McInnes 
176dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDestroy(PC);
177dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetFromOptions(PC);
178a313700dSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetType(PC,const PCType*);
17914c91fddSBarry Smith 
180a4fd02acSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatrix(PC,Mat*);
181dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
182dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
1835b116368SBarry Smith 
184dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetOperators(PC,Mat,Mat,MatStructure);
185dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperators(PC,Mat*,Mat*,MatStructure*);
186906ed7ccSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOperatorsSet(PC,PetscTruth*,PetscTruth*);
1874b0e389bSBarry Smith 
188dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCView(PC,PetscViewer);
1897bc3d0afSSatish Balay 
190dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetOptionsPrefix(PC,const char[]);
191dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCAppendOptionsPrefix(PC,const char[]);
1922dcb1b2aSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetOptionsPrefix(PC,const char*[]);
1938ed539a5SBarry Smith 
194dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCComputeExplicitOperator(PC,Mat*);
19571601f6fSBarry Smith 
196d6913704SBarry Smith /*
197d6913704SBarry Smith       These are used to provide extra scaling of preconditioned
1980f3b3ca1SHong Zhang    operator for time-stepping schemes like in SUNDIALS
199d6913704SBarry Smith */
200dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScale(PC,PetscTruth*);
201dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleLeft(PC,Vec,Vec);
202dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleRight(PC,Vec,Vec);
203dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCDiagonalScaleSet(PC,Vec);
204d6913704SBarry Smith 
20584cb2905SBarry Smith /* ------------- options specific to particular preconditioners --------- */
206329f5518SBarry Smith 
207dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowMax(PC);
208d87ca0b3SLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseRowSum(PC);
209cd47f5d9SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCJacobiSetUseAbs(PC);
210dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric(PC,MatSORType);
211dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega(PC,PetscReal);
212dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations(PC,PetscInt,PetscInt);
213d03aef70SBarry Smith 
214dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatSetOmega(PC,PetscReal);
215dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCEisenstatNoDiagonalScaling(PC);
216421c37bdSBarry Smith 
21715aa81f8SBarry Smith #define USE_PRECONDITIONER_MATRIX 0
21815aa81f8SBarry Smith #define USE_TRUE_MATRIX           1
219dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetUseTrueLocal(PC);
220dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
221dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
2221eb62cbbSBarry Smith 
223dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPSetUseTrue(PC);
224981c4779SBarry Smith 
2256891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
2266891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
2276891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
2286891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
2297319c654SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscTruth,PetscInt*,PCRichardsonConvergedReason*));
2306891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
2316891c3e4SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
232be29d3c6SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetContext(PC,void**);
233be29d3c6SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetContext(PC,void*);
234dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetName(PC,const char[]);
235dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellGetName(PC,char*[]);
236aabeff55SBarry Smith 
237dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot(PC,PetscReal);
238d90ac83dSHong Zhang 
239d90ac83dSHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftType(PC,MatFactorShiftType);
240d90ac83dSHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftAmount(PC,PetscReal);
241d90ac83dSHong Zhang 
242c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
24335bd34faSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
2442401956bSBarry Smith 
24555ba2a51SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC,PetscReal);
2468ff23777SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetColumnPivot(PC,PetscReal);
2472401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
248421c37bdSBarry Smith 
249e36ff184SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC,const MatOrderingType);
2502401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC,PetscTruth);
2512401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC,PetscTruth);
2522401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC);
25312a54f56SSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC);
2542401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC,PetscTruth);
255f5a88c2aSLois Curfman McInnes 
2562401956bSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC,PetscInt);
257b7c853c4SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
258b35a507dSBarry Smith 
2592b691e39SMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
2602b691e39SMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
261dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetOverlap(PC,PetscInt);
2626ed231c7SMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetSortIndices(PC,PetscTruth);
2633d957683SBarry Smith /*E
2643d957683SBarry Smith     PCASMType - Type of additive Schwarz method to use
2653d957683SBarry Smith 
2663d957683SBarry Smith $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
2673d957683SBarry Smith $                 and computed values in ghost regions are added together. Classical
2683d957683SBarry Smith $                 standard additive Schwarz
2693d957683SBarry Smith $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
2703d957683SBarry Smith $                    region are discarded. Default
2713d957683SBarry Smith $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
2723d957683SBarry Smith $                       region are added back in
2733d957683SBarry Smith $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
2743d957683SBarry Smith $                not very good.
2753d957683SBarry Smith 
2763d957683SBarry Smith    Level: beginner
2773d957683SBarry Smith 
2783d957683SBarry Smith .seealso: PCASMSetType()
2793d957683SBarry Smith E*/
280d252947aSBarry Smith typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
2819dcbbd2bSBarry Smith extern const char *PCASMTypes[];
2829dcbbd2bSBarry Smith 
283dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetType(PC,PCASMType);
28492bb6962SLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
2852b691e39SMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMDestroySubdomains(PetscInt,IS[],IS[]);
286dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt *,IS **);
287dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetUseInPlace(PC);
2882b691e39SMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
289dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
290981c4779SBarry Smith 
2913d957683SBarry Smith /*E
2923d957683SBarry Smith     PCCompositeType - Determines how two or more preconditioner are composed
2933d957683SBarry Smith 
2943d957683SBarry Smith $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
2953d957683SBarry Smith $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
2963d957683SBarry Smith $                                computed after the previous preconditioner application
297ccb205f8SBarry Smith $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
298ccb205f8SBarry Smith $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
2993d957683SBarry Smith $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
3003d957683SBarry Smith $                         where first preconditioner is built from alpha I + S and second from
3013d957683SBarry Smith $                         alpha I + R
3023d957683SBarry Smith 
3033d957683SBarry Smith    Level: beginner
3043d957683SBarry Smith 
3053d957683SBarry Smith .seealso: PCCompositeSetType()
3063d957683SBarry Smith E*/
3073b224e63SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
3089dcbbd2bSBarry Smith extern const char *PCCompositeTypes[];
3099dcbbd2bSBarry Smith 
310dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC);
311dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC,PCCompositeType);
312dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC,PCType);
3139a05e725SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC,PetscInt,PC *);
314dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC,PetscScalar);
315981c4779SBarry Smith 
31609a6bc64SHong Zhang EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber(PC,PetscInt);
317dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter(PC,VecScatter,VecScatter);
318dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators(PC,Mat*,Mat*);
319dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC(PC,PC*);
320da3a660dSBarry Smith 
321dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetEpsilon(PC,double);
322dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetNBSteps(PC,PetscInt);
323dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMax(PC,PetscInt);
324dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMaxNew(PC,PetscInt);
325dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetBlockSize(PC,PetscInt);
326dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetCacheSize(PC,PetscInt);
327dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetVerbose(PC,PetscInt);
328dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetSp(PC,PetscInt);
3293304466cSBarry Smith 
330dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC,const char[]);
331dd896f7dSSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC,const char*[]);
332dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
333dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
3343304466cSBarry Smith 
335dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC,PetscInt,PetscInt*);
336dba47a55SKris Buschelman EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC,PCCompositeType);
33751f519a2SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC,PetscInt);
33890fdf44cSMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC,IS);
339084e4875SJed Brown 
340084e4875SJed Brown /*E
341084e4875SJed Brown     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
342084e4875SJed Brown 
343084e4875SJed Brown     Level: intermediate
344084e4875SJed Brown 
345084e4875SJed Brown .seealso: PCFieldSplitSchurPrecondition()
346084e4875SJed Brown E*/
347084e4875SJed Brown typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_DIAG,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
348084e4875SJed Brown extern const char *PCFieldSplitSchurPreTypes[];
349084e4875SJed Brown 
350084e4875SJed Brown EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
35130ad9308SMatthew Knepley EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
3523d30b1ffSBarry Smith 
3532a6744ebSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetRestriction(PC,Mat);
3542a6744ebSBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetInterpolation(PC,Mat);
3552a6744ebSBarry Smith 
3561a3fd43dSSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetCoordinates(PC,PetscInt,PetscReal*);
3572102ba4dSSatish Balay EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSASetVectors(PC,PetscInt,PetscReal *);
3582102ba4dSSatish Balay 
3591d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPythonSetType(PC,const char[]);
36067fac13cSBarry Smith 
3616c699258SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCSetDM(PC,DM);
3626c699258SBarry Smith EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCGetDM(PC,DM*);
3636c699258SBarry Smith 
3646c699258SBarry Smith 
365e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
366d3b1e0e7SMatthew Knepley 
367123ea438SMatthew Knepley #endif /* __PETSCPC_H */
368