xref: /petsc/include/petscpc.h (revision 37f802242ef9160e6bb9371c3c5cccfcf23701e6)
1d03aef70SBarry Smith /*
237f753daSBarry Smith       Preconditioner module.
3d03aef70SBarry Smith */
40a835dfdSSatish Balay #if !defined(__PETSCPC_H)
50a835dfdSSatish Balay #define __PETSCPC_H
6e1589f56SBarry Smith #include "petscdm.h"
7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
8d03aef70SBarry Smith 
97087cfbeSBarry Smith extern PetscErrorCode   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*/
39a313700dSBarry Smith #define PCType 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"
51ab3e923fSDmitry Karpeev #define PCGASM            "gasm"
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"
61*37f80224SJose E. Roman #define PCPARMS           "parms"
620971522cSBarry Smith #define PCFIELDSPLIT      "fieldsplit"
63be16f70fSSatish Balay #define PCTFS             "tfs"
645582bec1SHong Zhang #define PCML              "ml"
6536a49750SBarry Smith #define PCPROMETHEUS      "prometheus"
662a6744ebSBarry Smith #define PCGALERKIN        "galerkin"
677233f9f0SBarry Smith #define PCEXOTIC          "exotic"
68161c5ca9SBarry Smith #define PCOPENMP          "openmp"
69cf037197Ssdaitch #define PCSUPPORTGRAPH    "supportgraph"
70f4b8409dSBarry Smith #define PCASA             "asa"
7124c02a0fSBarry Smith #define PCCP              "cp"
72628b8437SMatthew Knepley #define PCBFBT            "bfbt"
73519d70e2SJed Brown #define PCLSC             "lsc"
741d6018f0SLisandro Dalcin #define PCPYTHON          "python"
75f91d8e95SBarry Smith #define PCPFMG            "pfmg"
76d851a50bSGlenn Hammond #define PCSYSPFMG         "syspfmg"
77df826632SBarry Smith #define PCREDISTRIBUTE    "redistribute"
788154be41SBarry Smith #define PCSACUSP          "sacusp"
79bfc29b71SVictor Minden #define PCSACUSPPOLY      "sacusppoly"
808154be41SBarry Smith #define PCBICGSTABCUSP    "bicgstabcusp"
8127c67122SBarry Smith #define PCSVD             "svd"
8204b59e88SVictor Minden #define PCAINVCUSP        "ainvcusp"
83123ea438SMatthew Knepley 
84123ea438SMatthew Knepley /* Logging support */
857087cfbeSBarry Smith extern PetscClassId  PC_CLASSID;
86123ea438SMatthew Knepley 
873d957683SBarry Smith /*E
883d957683SBarry Smith     PCSide - If the preconditioner is to be applied to the left, right
893d957683SBarry Smith      or symmetrically around the operator.
90d03aef70SBarry Smith 
913d957683SBarry Smith    Level: beginner
9221e95762SBarry Smith 
933d957683SBarry Smith .seealso:
943d957683SBarry Smith E*/
95aad2872bSLois Curfman McInnes typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
969dcbbd2bSBarry Smith extern const char *PCSides[];
9772b7852fSLois Curfman McInnes 
987087cfbeSBarry Smith extern PetscErrorCode  PCCreate(MPI_Comm,PC*);
997087cfbeSBarry Smith extern PetscErrorCode  PCSetType(PC,const PCType);
1007087cfbeSBarry Smith extern PetscErrorCode  PCSetUp(PC);
1017087cfbeSBarry Smith extern PetscErrorCode  PCSetUpOnBlocks(PC);
1027087cfbeSBarry Smith extern PetscErrorCode  PCApply(PC,Vec,Vec);
1037087cfbeSBarry Smith extern PetscErrorCode  PCApplySymmetricLeft(PC,Vec,Vec);
1047087cfbeSBarry Smith extern PetscErrorCode  PCApplySymmetricRight(PC,Vec,Vec);
1057087cfbeSBarry Smith extern PetscErrorCode  PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
1067087cfbeSBarry Smith extern PetscErrorCode  PCApplyTranspose(PC,Vec,Vec);
1077087cfbeSBarry Smith extern PetscErrorCode  PCApplyTransposeExists(PC,PetscBool *);
1087087cfbeSBarry Smith extern PetscErrorCode  PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
1094d0a8057SBarry Smith 
1104d0a8057SBarry Smith /*E
1114d0a8057SBarry Smith     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
1124d0a8057SBarry Smith 
1134d0a8057SBarry Smith    Level: advanced
1144d0a8057SBarry Smith 
1154d0a8057SBarry Smith    Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
1164d0a8057SBarry Smith 
1174d0a8057SBarry Smith .seealso: PCApplyRichardson()
1184d0a8057SBarry Smith E*/
1194d0a8057SBarry Smith typedef enum {
1204d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_RTOL               =  2,
1214d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_ATOL               =  3,
1224d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_ITS                =  4,
1234d0a8057SBarry Smith               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
1244d0a8057SBarry Smith 
1257087cfbeSBarry Smith extern PetscErrorCode  PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
1267087cfbeSBarry Smith extern PetscErrorCode  PCApplyRichardsonExists(PC,PetscBool *);
1277087cfbeSBarry Smith extern PetscErrorCode  PCSetInitialGuessNonzero(PC,PetscBool );
12884cb2905SBarry Smith 
1297087cfbeSBarry Smith extern PetscErrorCode  PCRegisterDestroy(void);
1307087cfbeSBarry Smith extern PetscErrorCode  PCRegisterAll(const char[]);
131ace3abfcSBarry Smith extern PetscBool  PCRegisterAllCalled;
13284cb2905SBarry Smith 
1337087cfbeSBarry Smith extern PetscErrorCode  PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
13430de9b25SBarry Smith 
13530de9b25SBarry Smith /*MC
13630de9b25SBarry Smith    PCRegisterDynamic - Adds a method to the preconditioner package.
13730de9b25SBarry Smith 
13830de9b25SBarry Smith    Synopsis:
1391890ba74SBarry Smith    PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC))
14030de9b25SBarry Smith 
14130de9b25SBarry Smith    Not collective
14230de9b25SBarry Smith 
14330de9b25SBarry Smith    Input Parameters:
14430de9b25SBarry Smith +  name_solver - name of a new user-defined solver
14530de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
14630de9b25SBarry Smith .  name_create - name of routine to create method context
14730de9b25SBarry Smith -  routine_create - routine to create method context
14830de9b25SBarry Smith 
14930de9b25SBarry Smith    Notes:
15030de9b25SBarry Smith    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
15130de9b25SBarry Smith 
15230de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
15330de9b25SBarry Smith    is ignored.
15430de9b25SBarry Smith 
15530de9b25SBarry Smith    Sample usage:
15630de9b25SBarry Smith .vb
15730de9b25SBarry Smith    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
15830de9b25SBarry Smith               "MySolverCreate",MySolverCreate);
15930de9b25SBarry Smith .ve
16030de9b25SBarry Smith 
16130de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
16230de9b25SBarry Smith $     PCSetType(pc,"my_solver")
16330de9b25SBarry Smith    or at runtime via the option
16430de9b25SBarry Smith $     -pc_type my_solver
16530de9b25SBarry Smith 
16630de9b25SBarry Smith    Level: advanced
16730de9b25SBarry Smith 
168ab901514SBarry Smith    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
16930de9b25SBarry Smith            occuring in pathname will be replaced with appropriate values.
17030de9b25SBarry Smith          If your function is not being put into a shared library then use PCRegister() instead
17130de9b25SBarry Smith 
17230de9b25SBarry Smith .keywords: PC, register
17330de9b25SBarry Smith 
17430de9b25SBarry Smith .seealso: PCRegisterAll(), PCRegisterDestroy()
17530de9b25SBarry Smith M*/
176aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
177f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
1786df38c32SLois Curfman McInnes #else
179f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
1806df38c32SLois Curfman McInnes #endif
1816df38c32SLois Curfman McInnes 
182cd8bf0efSSatish Balay extern PetscErrorCode  PCReset(PC);
1837087cfbeSBarry Smith extern PetscErrorCode  PCDestroy(PC);
1847087cfbeSBarry Smith extern PetscErrorCode  PCSetFromOptions(PC);
1857087cfbeSBarry Smith extern PetscErrorCode  PCGetType(PC,const PCType*);
18614c91fddSBarry Smith 
1877087cfbeSBarry Smith extern PetscErrorCode  PCFactorGetMatrix(PC,Mat*);
1887087cfbeSBarry Smith extern PetscErrorCode  PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
1897087cfbeSBarry Smith extern PetscErrorCode  PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
1905b116368SBarry Smith 
1917087cfbeSBarry Smith extern PetscErrorCode  PCSetOperators(PC,Mat,Mat,MatStructure);
1927087cfbeSBarry Smith extern PetscErrorCode  PCGetOperators(PC,Mat*,Mat*,MatStructure*);
1937087cfbeSBarry Smith extern PetscErrorCode  PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
1944b0e389bSBarry Smith 
1957087cfbeSBarry Smith extern PetscErrorCode  PCView(PC,PetscViewer);
1967bc3d0afSSatish Balay 
1977087cfbeSBarry Smith extern PetscErrorCode  PCSetOptionsPrefix(PC,const char[]);
1987087cfbeSBarry Smith extern PetscErrorCode  PCAppendOptionsPrefix(PC,const char[]);
1997087cfbeSBarry Smith extern PetscErrorCode  PCGetOptionsPrefix(PC,const char*[]);
2008ed539a5SBarry Smith 
2017087cfbeSBarry Smith extern PetscErrorCode  PCComputeExplicitOperator(PC,Mat*);
20271601f6fSBarry Smith 
203d6913704SBarry Smith /*
204d6913704SBarry Smith       These are used to provide extra scaling of preconditioned
2050f3b3ca1SHong Zhang    operator for time-stepping schemes like in SUNDIALS
206d6913704SBarry Smith */
2077087cfbeSBarry Smith extern PetscErrorCode  PCGetDiagonalScale(PC,PetscBool *);
2087087cfbeSBarry Smith extern PetscErrorCode  PCDiagonalScaleLeft(PC,Vec,Vec);
2097087cfbeSBarry Smith extern PetscErrorCode  PCDiagonalScaleRight(PC,Vec,Vec);
2107087cfbeSBarry Smith extern PetscErrorCode  PCSetDiagonalScale(PC,Vec);
211d6913704SBarry Smith 
21284cb2905SBarry Smith /* ------------- options specific to particular preconditioners --------- */
213329f5518SBarry Smith 
2147087cfbeSBarry Smith extern PetscErrorCode  PCJacobiSetUseRowMax(PC);
2157087cfbeSBarry Smith extern PetscErrorCode  PCJacobiSetUseRowSum(PC);
2167087cfbeSBarry Smith extern PetscErrorCode  PCJacobiSetUseAbs(PC);
2177087cfbeSBarry Smith extern PetscErrorCode  PCSORSetSymmetric(PC,MatSORType);
2187087cfbeSBarry Smith extern PetscErrorCode  PCSORSetOmega(PC,PetscReal);
2197087cfbeSBarry Smith extern PetscErrorCode  PCSORSetIterations(PC,PetscInt,PetscInt);
220d03aef70SBarry Smith 
2217087cfbeSBarry Smith extern PetscErrorCode  PCEisenstatSetOmega(PC,PetscReal);
2227087cfbeSBarry Smith extern PetscErrorCode  PCEisenstatNoDiagonalScaling(PC);
223421c37bdSBarry Smith 
22415aa81f8SBarry Smith #define USE_PRECONDITIONER_MATRIX 0
22515aa81f8SBarry Smith #define USE_TRUE_MATRIX           1
2267087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiSetUseTrueLocal(PC);
2277087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
2287087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
2291eb62cbbSBarry Smith 
2307087cfbeSBarry Smith extern PetscErrorCode  PCKSPSetUseTrue(PC);
231981c4779SBarry Smith 
2327087cfbeSBarry Smith extern PetscErrorCode  PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
2337087cfbeSBarry Smith extern PetscErrorCode  PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
2347087cfbeSBarry Smith extern PetscErrorCode  PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
2357087cfbeSBarry Smith extern PetscErrorCode  PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
2367087cfbeSBarry Smith extern PetscErrorCode  PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
2377087cfbeSBarry Smith extern PetscErrorCode  PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
2387087cfbeSBarry Smith extern PetscErrorCode  PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
2397087cfbeSBarry Smith extern PetscErrorCode  PCShellGetContext(PC,void**);
2407087cfbeSBarry Smith extern PetscErrorCode  PCShellSetContext(PC,void*);
2417087cfbeSBarry Smith extern PetscErrorCode  PCShellSetName(PC,const char[]);
2427087cfbeSBarry Smith extern PetscErrorCode  PCShellGetName(PC,char*[]);
243aabeff55SBarry Smith 
2447087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetZeroPivot(PC,PetscReal);
245d90ac83dSHong Zhang 
2467087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetShiftType(PC,MatFactorShiftType);
2477087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetShiftAmount(PC,PetscReal);
248d90ac83dSHong Zhang 
2497087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
2507087cfbeSBarry Smith extern PetscErrorCode  PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
251f8260c8fSBarry Smith extern PetscErrorCode  PCFactorSetUpMatSolverPackage(PC);
2522401956bSBarry Smith 
2537087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetFill(PC,PetscReal);
2547087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetColumnPivot(PC,PetscReal);
2557087cfbeSBarry Smith extern PetscErrorCode  PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
256421c37bdSBarry Smith 
2577087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetMatOrderingType(PC,const MatOrderingType);
2587087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetReuseOrdering(PC,PetscBool );
2597087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetReuseFill(PC,PetscBool );
2607087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetUseInPlace(PC);
2617087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetAllowDiagonalFill(PC);
2627087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetPivotInBlocks(PC,PetscBool );
263f5a88c2aSLois Curfman McInnes 
2647087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetLevels(PC,PetscInt);
2657087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
266b35a507dSBarry Smith 
2677087cfbeSBarry Smith extern PetscErrorCode  PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
2687087cfbeSBarry Smith extern PetscErrorCode  PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
2697087cfbeSBarry Smith extern PetscErrorCode  PCASMSetOverlap(PC,PetscInt);
2707087cfbeSBarry Smith extern PetscErrorCode  PCASMSetSortIndices(PC,PetscBool );
271f746d493SDmitry Karpeev 
2723d957683SBarry Smith /*E
2733d957683SBarry Smith     PCASMType - Type of additive Schwarz method to use
2743d957683SBarry Smith 
2753d957683SBarry Smith $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
2763d957683SBarry Smith $                 and computed values in ghost regions are added together. Classical
2773d957683SBarry Smith $                 standard additive Schwarz
2783d957683SBarry Smith $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
2793d957683SBarry Smith $                    region are discarded. Default
2803d957683SBarry Smith $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
2813d957683SBarry Smith $                       region are added back in
2823d957683SBarry Smith $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
2833d957683SBarry Smith $                not very good.
2843d957683SBarry Smith 
2853d957683SBarry Smith    Level: beginner
2863d957683SBarry Smith 
2873d957683SBarry Smith .seealso: PCASMSetType()
2883d957683SBarry Smith E*/
289d252947aSBarry Smith typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
2909dcbbd2bSBarry Smith extern const char *PCASMTypes[];
2919dcbbd2bSBarry Smith 
2927087cfbeSBarry Smith extern PetscErrorCode  PCASMSetType(PC,PCASMType);
2937087cfbeSBarry Smith extern PetscErrorCode  PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
2947087cfbeSBarry Smith extern PetscErrorCode  PCASMDestroySubdomains(PetscInt,IS[],IS[]);
2957087cfbeSBarry Smith extern PetscErrorCode  PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
2967087cfbeSBarry Smith extern PetscErrorCode  PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
2977087cfbeSBarry Smith extern PetscErrorCode  PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
298981c4779SBarry Smith 
2993d957683SBarry Smith /*E
300f746d493SDmitry Karpeev     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per domain)
301f746d493SDmitry Karpeev 
302f746d493SDmitry Karpeev $  PC_GASM_BASIC - symmetric version where residuals from the ghost points are used
303f746d493SDmitry Karpeev $                 and computed values in ghost regions are added together. Classical
304f746d493SDmitry Karpeev $                 standard additive Schwarz
305f746d493SDmitry Karpeev $  PC_GASM_RESTRICT - residuals from ghost points are used but computed values in ghost
306f746d493SDmitry Karpeev $                    region are discarded. Default
307f746d493SDmitry Karpeev $  PC_GASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
308f746d493SDmitry Karpeev $                       region are added back in
309f746d493SDmitry Karpeev $  PC_GASM_NONE - ghost point residuals are not used, computed ghost values are discarded
310f746d493SDmitry Karpeev $                not very good.
311f746d493SDmitry Karpeev 
312f746d493SDmitry Karpeev    Level: beginner
313f746d493SDmitry Karpeev 
314f746d493SDmitry Karpeev .seealso: PCGASMSetType()
315f746d493SDmitry Karpeev E*/
316f746d493SDmitry Karpeev typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
317f746d493SDmitry Karpeev extern const char *PCGASMTypes[];
318f746d493SDmitry Karpeev 
3197087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
3207087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetTotalSubdomains(PC,PetscInt);
3217087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetOverlap(PC,PetscInt);
3227087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetSortIndices(PC,PetscBool );
323f746d493SDmitry Karpeev 
3247087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetType(PC,PCGASMType);
3257087cfbeSBarry Smith extern PetscErrorCode  PCGASMCreateSubdomains(Mat,PetscInt,IS*[]);
3267087cfbeSBarry Smith extern PetscErrorCode  PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
3277087cfbeSBarry Smith extern PetscErrorCode  PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
3287087cfbeSBarry Smith extern PetscErrorCode  PCGASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
3297087cfbeSBarry Smith extern PetscErrorCode  PCGASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
330f746d493SDmitry Karpeev 
331f746d493SDmitry Karpeev /*E
3323d957683SBarry Smith     PCCompositeType - Determines how two or more preconditioner are composed
3333d957683SBarry Smith 
3343d957683SBarry Smith $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
3353d957683SBarry Smith $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
3363d957683SBarry Smith $                                computed after the previous preconditioner application
337ccb205f8SBarry Smith $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
338ccb205f8SBarry Smith $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
3393d957683SBarry Smith $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
3403d957683SBarry Smith $                         where first preconditioner is built from alpha I + S and second from
3413d957683SBarry Smith $                         alpha I + R
3423d957683SBarry Smith 
3433d957683SBarry Smith    Level: beginner
3443d957683SBarry Smith 
3453d957683SBarry Smith .seealso: PCCompositeSetType()
3463d957683SBarry Smith E*/
3473b224e63SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
3489dcbbd2bSBarry Smith extern const char *PCCompositeTypes[];
3499dcbbd2bSBarry Smith 
3507087cfbeSBarry Smith extern PetscErrorCode  PCCompositeSetUseTrue(PC);
3517087cfbeSBarry Smith extern PetscErrorCode  PCCompositeSetType(PC,PCCompositeType);
3527087cfbeSBarry Smith extern PetscErrorCode  PCCompositeAddPC(PC,PCType);
3537087cfbeSBarry Smith extern PetscErrorCode  PCCompositeGetPC(PC,PetscInt,PC *);
3547087cfbeSBarry Smith extern PetscErrorCode  PCCompositeSpecialSetAlpha(PC,PetscScalar);
355981c4779SBarry Smith 
3567087cfbeSBarry Smith extern PetscErrorCode  PCRedundantSetNumber(PC,PetscInt);
3577087cfbeSBarry Smith extern PetscErrorCode  PCRedundantSetScatter(PC,VecScatter,VecScatter);
3587087cfbeSBarry Smith extern PetscErrorCode  PCRedundantGetOperators(PC,Mat*,Mat*);
359da3a660dSBarry Smith 
3607087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetEpsilon(PC,double);
3617087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetNBSteps(PC,PetscInt);
3627087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetMax(PC,PetscInt);
3637087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetMaxNew(PC,PetscInt);
3647087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetBlockSize(PC,PetscInt);
3657087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetCacheSize(PC,PetscInt);
3667087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetVerbose(PC,PetscInt);
3677087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetSp(PC,PetscInt);
3683304466cSBarry Smith 
3697087cfbeSBarry Smith extern PetscErrorCode  PCHYPRESetType(PC,const char[]);
3707087cfbeSBarry Smith extern PetscErrorCode  PCHYPREGetType(PC,const char*[]);
3717087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
3727087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
3733304466cSBarry Smith 
3747087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*);
3757087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSetType(PC,PCCompositeType);
3767087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSetBlockSize(PC,PetscInt);
3777087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSetIS(PC,const char[],IS);
37857a9adfeSMatthew G Knepley extern PetscErrorCode  PCFieldSplitGetIS(PC,const char[],IS*);
379084e4875SJed Brown 
380084e4875SJed Brown /*E
381084e4875SJed Brown     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
382084e4875SJed Brown 
383084e4875SJed Brown     Level: intermediate
384084e4875SJed Brown 
385084e4875SJed Brown .seealso: PCFieldSplitSchurPrecondition()
386084e4875SJed Brown E*/
387084e4875SJed Brown typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_DIAG,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
388f7c28744SJed Brown extern const char *const PCFieldSplitSchurPreTypes[];
389084e4875SJed Brown 
3907087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
3917087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
3923d30b1ffSBarry Smith 
3937087cfbeSBarry Smith extern PetscErrorCode  PCGalerkinSetRestriction(PC,Mat);
3947087cfbeSBarry Smith extern PetscErrorCode  PCGalerkinSetInterpolation(PC,Mat);
3952a6744ebSBarry Smith 
3967087cfbeSBarry Smith extern PetscErrorCode  PCSetCoordinates(PC,PetscInt,PetscReal*);
3977087cfbeSBarry Smith extern PetscErrorCode  PCSASetVectors(PC,PetscInt,PetscReal *);
3982102ba4dSSatish Balay 
3997087cfbeSBarry Smith extern PetscErrorCode  PCPythonSetType(PC,const char[]);
40067fac13cSBarry Smith 
4017087cfbeSBarry Smith extern PetscErrorCode  PCSetDM(PC,DM);
4027087cfbeSBarry Smith extern PetscErrorCode  PCGetDM(PC,DM*);
4036c699258SBarry Smith 
4048154be41SBarry Smith extern PetscErrorCode  PCBiCGStabCUSPSetTolerance(PC,PetscReal);
4053d740270SVictor Minden extern PetscErrorCode  PCBiCGStabCUSPSetIterations(PC,PetscInt);
4063d740270SVictor Minden extern PetscErrorCode  PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
4076c699258SBarry Smith 
4086de8952cSVictor Minden extern PetscErrorCode  PCAINVCUSPSetDropTolerance(PC,PetscReal);
4096de8952cSVictor Minden extern PetscErrorCode  PCAINVCUSPUseScaling(PC,PetscBool);
4106de8952cSVictor Minden extern PetscErrorCode  PCAINVCUSPSetNonzeros(PC,PetscInt);
411adbda850SVictor Minden extern PetscErrorCode  PCAINVCUSPSetLinParameter(PC,PetscInt);
412*37f80224SJose E. Roman /*E
413*37f80224SJose E. Roman     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
414*37f80224SJose E. Roman 
415*37f80224SJose E. Roman     Level: intermediate
416*37f80224SJose E. Roman 
417*37f80224SJose E. Roman .seealso: PCPARMSSetGlobal()
418*37f80224SJose E. Roman E*/
419*37f80224SJose E. Roman typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
420*37f80224SJose E. Roman extern const char *PCPARMSGlobalTypes[];
421*37f80224SJose E. Roman /*E
422*37f80224SJose E. Roman     PCPARMSLocalType - Determines the local preconditioner method in PARMS
423*37f80224SJose E. Roman 
424*37f80224SJose E. Roman     Level: intermediate
425*37f80224SJose E. Roman 
426*37f80224SJose E. Roman .seealso: PCPARMSSetLocal()
427*37f80224SJose E. Roman E*/
428*37f80224SJose E. Roman typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
429*37f80224SJose E. Roman extern const char *PCPARMSLocalTypes[];
430*37f80224SJose E. Roman 
431*37f80224SJose E. Roman extern PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type);
432*37f80224SJose E. Roman extern PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type);
433*37f80224SJose E. Roman extern PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits);
434*37f80224SJose E. Roman extern PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart);
435*37f80224SJose E. Roman extern PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym);
436*37f80224SJose E. Roman extern PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2);
437*37f80224SJose E. Roman 
438e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
439d3b1e0e7SMatthew Knepley 
440123ea438SMatthew Knepley #endif /* __PETSCPC_H */
441