xref: /petsc/include/petscpc.h (revision cd8bf0ef4edbe2f9530b1f349b342a04731c463c)
1d03aef70SBarry Smith /*
237f753daSBarry Smith       Preconditioner module.
3d03aef70SBarry Smith */
40a835dfdSSatish Balay #if !defined(__PETSCPC_H)
50a835dfdSSatish Balay #define __PETSCPC_H
60a835dfdSSatish Balay #include "petscmat.h"
7e1589f56SBarry Smith #include "petscdm.h"
8e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
9d03aef70SBarry Smith 
107087cfbeSBarry Smith extern PetscErrorCode   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"
52ab3e923fSDmitry Karpeev #define PCGASM            "gasm"
5394b7f48cSBarry Smith #define PCKSP             "ksp"
5482bf6240SBarry Smith #define PCCOMPOSITE       "composite"
55421c37bdSBarry Smith #define PCREDUNDANT       "redundant"
5627b520f0SBarry Smith #define PCSPAI            "spai"
57186905e3SBarry Smith #define PCNN              "nn"
584bbc92c1SBarry Smith #define PCCHOLESKY        "cholesky"
593a7fca6bSBarry Smith #define PCPBJACOBI        "pbjacobi"
607f5ff6fdSBarry Smith #define PCMAT             "mat"
61c4888f26SBarry Smith #define PCHYPRE           "hypre"
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"
798154be41SBarry Smith #define PCBICGSTABCUSP    "bicgstabcusp"
8027c67122SBarry Smith #define PCSVD             "svd"
81123ea438SMatthew Knepley 
82123ea438SMatthew Knepley /* Logging support */
837087cfbeSBarry Smith extern PetscClassId  PC_CLASSID;
84123ea438SMatthew Knepley 
853d957683SBarry Smith /*E
863d957683SBarry Smith     PCSide - If the preconditioner is to be applied to the left, right
873d957683SBarry Smith      or symmetrically around the operator.
88d03aef70SBarry Smith 
893d957683SBarry Smith    Level: beginner
9021e95762SBarry Smith 
913d957683SBarry Smith .seealso:
923d957683SBarry Smith E*/
93aad2872bSLois Curfman McInnes typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
949dcbbd2bSBarry Smith extern const char *PCSides[];
9572b7852fSLois Curfman McInnes 
967087cfbeSBarry Smith extern PetscErrorCode  PCCreate(MPI_Comm,PC*);
977087cfbeSBarry Smith extern PetscErrorCode  PCSetType(PC,const PCType);
987087cfbeSBarry Smith extern PetscErrorCode  PCSetUp(PC);
997087cfbeSBarry Smith extern PetscErrorCode  PCSetUpOnBlocks(PC);
1007087cfbeSBarry Smith extern PetscErrorCode  PCApply(PC,Vec,Vec);
1017087cfbeSBarry Smith extern PetscErrorCode  PCApplySymmetricLeft(PC,Vec,Vec);
1027087cfbeSBarry Smith extern PetscErrorCode  PCApplySymmetricRight(PC,Vec,Vec);
1037087cfbeSBarry Smith extern PetscErrorCode  PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
1047087cfbeSBarry Smith extern PetscErrorCode  PCApplyTranspose(PC,Vec,Vec);
1057087cfbeSBarry Smith extern PetscErrorCode  PCApplyTransposeExists(PC,PetscBool *);
1067087cfbeSBarry Smith extern PetscErrorCode  PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
1074d0a8057SBarry Smith 
1084d0a8057SBarry Smith /*E
1094d0a8057SBarry Smith     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
1104d0a8057SBarry Smith 
1114d0a8057SBarry Smith    Level: advanced
1124d0a8057SBarry Smith 
1134d0a8057SBarry Smith    Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
1144d0a8057SBarry Smith 
1154d0a8057SBarry Smith .seealso: PCApplyRichardson()
1164d0a8057SBarry Smith E*/
1174d0a8057SBarry Smith typedef enum {
1184d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_RTOL               =  2,
1194d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_ATOL               =  3,
1204d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_ITS                =  4,
1214d0a8057SBarry Smith               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
1224d0a8057SBarry Smith 
1237087cfbeSBarry Smith extern PetscErrorCode  PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
1247087cfbeSBarry Smith extern PetscErrorCode  PCApplyRichardsonExists(PC,PetscBool *);
1257087cfbeSBarry Smith extern PetscErrorCode  PCSetInitialGuessNonzero(PC,PetscBool );
12684cb2905SBarry Smith 
1277087cfbeSBarry Smith extern PetscErrorCode  PCRegisterDestroy(void);
1287087cfbeSBarry Smith extern PetscErrorCode  PCRegisterAll(const char[]);
129ace3abfcSBarry Smith extern PetscBool  PCRegisterAllCalled;
13084cb2905SBarry Smith 
1317087cfbeSBarry Smith extern PetscErrorCode  PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
13230de9b25SBarry Smith 
13330de9b25SBarry Smith /*MC
13430de9b25SBarry Smith    PCRegisterDynamic - Adds a method to the preconditioner package.
13530de9b25SBarry Smith 
13630de9b25SBarry Smith    Synopsis:
1371890ba74SBarry Smith    PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC))
13830de9b25SBarry Smith 
13930de9b25SBarry Smith    Not collective
14030de9b25SBarry Smith 
14130de9b25SBarry Smith    Input Parameters:
14230de9b25SBarry Smith +  name_solver - name of a new user-defined solver
14330de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
14430de9b25SBarry Smith .  name_create - name of routine to create method context
14530de9b25SBarry Smith -  routine_create - routine to create method context
14630de9b25SBarry Smith 
14730de9b25SBarry Smith    Notes:
14830de9b25SBarry Smith    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
14930de9b25SBarry Smith 
15030de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
15130de9b25SBarry Smith    is ignored.
15230de9b25SBarry Smith 
15330de9b25SBarry Smith    Sample usage:
15430de9b25SBarry Smith .vb
15530de9b25SBarry Smith    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
15630de9b25SBarry Smith               "MySolverCreate",MySolverCreate);
15730de9b25SBarry Smith .ve
15830de9b25SBarry Smith 
15930de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
16030de9b25SBarry Smith $     PCSetType(pc,"my_solver")
16130de9b25SBarry Smith    or at runtime via the option
16230de9b25SBarry Smith $     -pc_type my_solver
16330de9b25SBarry Smith 
16430de9b25SBarry Smith    Level: advanced
16530de9b25SBarry Smith 
166ab901514SBarry Smith    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
16730de9b25SBarry Smith            occuring in pathname will be replaced with appropriate values.
16830de9b25SBarry Smith          If your function is not being put into a shared library then use PCRegister() instead
16930de9b25SBarry Smith 
17030de9b25SBarry Smith .keywords: PC, register
17130de9b25SBarry Smith 
17230de9b25SBarry Smith .seealso: PCRegisterAll(), PCRegisterDestroy()
17330de9b25SBarry Smith M*/
174aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
175f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
1766df38c32SLois Curfman McInnes #else
177f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
1786df38c32SLois Curfman McInnes #endif
1796df38c32SLois Curfman McInnes 
180*cd8bf0efSSatish Balay extern PetscErrorCode  PCReset(PC);
1817087cfbeSBarry Smith extern PetscErrorCode  PCDestroy(PC);
1827087cfbeSBarry Smith extern PetscErrorCode  PCSetFromOptions(PC);
1837087cfbeSBarry Smith extern PetscErrorCode  PCGetType(PC,const PCType*);
18414c91fddSBarry Smith 
1857087cfbeSBarry Smith extern PetscErrorCode  PCFactorGetMatrix(PC,Mat*);
1867087cfbeSBarry Smith extern PetscErrorCode  PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
1877087cfbeSBarry Smith extern PetscErrorCode  PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
1885b116368SBarry Smith 
1897087cfbeSBarry Smith extern PetscErrorCode  PCSetOperators(PC,Mat,Mat,MatStructure);
1907087cfbeSBarry Smith extern PetscErrorCode  PCGetOperators(PC,Mat*,Mat*,MatStructure*);
1917087cfbeSBarry Smith extern PetscErrorCode  PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
1924b0e389bSBarry Smith 
1937087cfbeSBarry Smith extern PetscErrorCode  PCView(PC,PetscViewer);
1947bc3d0afSSatish Balay 
1957087cfbeSBarry Smith extern PetscErrorCode  PCSetOptionsPrefix(PC,const char[]);
1967087cfbeSBarry Smith extern PetscErrorCode  PCAppendOptionsPrefix(PC,const char[]);
1977087cfbeSBarry Smith extern PetscErrorCode  PCGetOptionsPrefix(PC,const char*[]);
1988ed539a5SBarry Smith 
1997087cfbeSBarry Smith extern PetscErrorCode  PCComputeExplicitOperator(PC,Mat*);
20071601f6fSBarry Smith 
201d6913704SBarry Smith /*
202d6913704SBarry Smith       These are used to provide extra scaling of preconditioned
2030f3b3ca1SHong Zhang    operator for time-stepping schemes like in SUNDIALS
204d6913704SBarry Smith */
2057087cfbeSBarry Smith extern PetscErrorCode  PCGetDiagonalScale(PC,PetscBool *);
2067087cfbeSBarry Smith extern PetscErrorCode  PCDiagonalScaleLeft(PC,Vec,Vec);
2077087cfbeSBarry Smith extern PetscErrorCode  PCDiagonalScaleRight(PC,Vec,Vec);
2087087cfbeSBarry Smith extern PetscErrorCode  PCSetDiagonalScale(PC,Vec);
209d6913704SBarry Smith 
21084cb2905SBarry Smith /* ------------- options specific to particular preconditioners --------- */
211329f5518SBarry Smith 
2127087cfbeSBarry Smith extern PetscErrorCode  PCJacobiSetUseRowMax(PC);
2137087cfbeSBarry Smith extern PetscErrorCode  PCJacobiSetUseRowSum(PC);
2147087cfbeSBarry Smith extern PetscErrorCode  PCJacobiSetUseAbs(PC);
2157087cfbeSBarry Smith extern PetscErrorCode  PCSORSetSymmetric(PC,MatSORType);
2167087cfbeSBarry Smith extern PetscErrorCode  PCSORSetOmega(PC,PetscReal);
2177087cfbeSBarry Smith extern PetscErrorCode  PCSORSetIterations(PC,PetscInt,PetscInt);
218d03aef70SBarry Smith 
2197087cfbeSBarry Smith extern PetscErrorCode  PCEisenstatSetOmega(PC,PetscReal);
2207087cfbeSBarry Smith extern PetscErrorCode  PCEisenstatNoDiagonalScaling(PC);
221421c37bdSBarry Smith 
22215aa81f8SBarry Smith #define USE_PRECONDITIONER_MATRIX 0
22315aa81f8SBarry Smith #define USE_TRUE_MATRIX           1
2247087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiSetUseTrueLocal(PC);
2257087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
2267087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
2271eb62cbbSBarry Smith 
2287087cfbeSBarry Smith extern PetscErrorCode  PCKSPSetUseTrue(PC);
229981c4779SBarry Smith 
2307087cfbeSBarry Smith extern PetscErrorCode  PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
2317087cfbeSBarry Smith extern PetscErrorCode  PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
2327087cfbeSBarry Smith extern PetscErrorCode  PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
2337087cfbeSBarry Smith extern PetscErrorCode  PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
2347087cfbeSBarry Smith extern PetscErrorCode  PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
2357087cfbeSBarry Smith extern PetscErrorCode  PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
2367087cfbeSBarry Smith extern PetscErrorCode  PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
2377087cfbeSBarry Smith extern PetscErrorCode  PCShellGetContext(PC,void**);
2387087cfbeSBarry Smith extern PetscErrorCode  PCShellSetContext(PC,void*);
2397087cfbeSBarry Smith extern PetscErrorCode  PCShellSetName(PC,const char[]);
2407087cfbeSBarry Smith extern PetscErrorCode  PCShellGetName(PC,char*[]);
241aabeff55SBarry Smith 
2427087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetZeroPivot(PC,PetscReal);
243d90ac83dSHong Zhang 
2447087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetShiftType(PC,MatFactorShiftType);
2457087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetShiftAmount(PC,PetscReal);
246d90ac83dSHong Zhang 
2477087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
2487087cfbeSBarry Smith extern PetscErrorCode  PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
249f8260c8fSBarry Smith extern PetscErrorCode  PCFactorSetUpMatSolverPackage(PC);
2502401956bSBarry Smith 
2517087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetFill(PC,PetscReal);
2527087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetColumnPivot(PC,PetscReal);
2537087cfbeSBarry Smith extern PetscErrorCode  PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
254421c37bdSBarry Smith 
2557087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetMatOrderingType(PC,const MatOrderingType);
2567087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetReuseOrdering(PC,PetscBool );
2577087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetReuseFill(PC,PetscBool );
2587087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetUseInPlace(PC);
2597087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetAllowDiagonalFill(PC);
2607087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetPivotInBlocks(PC,PetscBool );
261f5a88c2aSLois Curfman McInnes 
2627087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetLevels(PC,PetscInt);
2637087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
264b35a507dSBarry Smith 
2657087cfbeSBarry Smith extern PetscErrorCode  PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
2667087cfbeSBarry Smith extern PetscErrorCode  PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
2677087cfbeSBarry Smith extern PetscErrorCode  PCASMSetOverlap(PC,PetscInt);
2687087cfbeSBarry Smith extern PetscErrorCode  PCASMSetSortIndices(PC,PetscBool );
269f746d493SDmitry Karpeev 
2703d957683SBarry Smith /*E
2713d957683SBarry Smith     PCASMType - Type of additive Schwarz method to use
2723d957683SBarry Smith 
2733d957683SBarry Smith $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
2743d957683SBarry Smith $                 and computed values in ghost regions are added together. Classical
2753d957683SBarry Smith $                 standard additive Schwarz
2763d957683SBarry Smith $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
2773d957683SBarry Smith $                    region are discarded. Default
2783d957683SBarry Smith $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
2793d957683SBarry Smith $                       region are added back in
2803d957683SBarry Smith $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
2813d957683SBarry Smith $                not very good.
2823d957683SBarry Smith 
2833d957683SBarry Smith    Level: beginner
2843d957683SBarry Smith 
2853d957683SBarry Smith .seealso: PCASMSetType()
2863d957683SBarry Smith E*/
287d252947aSBarry Smith typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
2889dcbbd2bSBarry Smith extern const char *PCASMTypes[];
2899dcbbd2bSBarry Smith 
2907087cfbeSBarry Smith extern PetscErrorCode  PCASMSetType(PC,PCASMType);
2917087cfbeSBarry Smith extern PetscErrorCode  PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
2927087cfbeSBarry Smith extern PetscErrorCode  PCASMDestroySubdomains(PetscInt,IS[],IS[]);
2937087cfbeSBarry Smith extern PetscErrorCode  PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
2947087cfbeSBarry Smith extern PetscErrorCode  PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
2957087cfbeSBarry Smith extern PetscErrorCode  PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
296981c4779SBarry Smith 
2973d957683SBarry Smith /*E
298f746d493SDmitry Karpeev     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per domain)
299f746d493SDmitry Karpeev 
300f746d493SDmitry Karpeev $  PC_GASM_BASIC - symmetric version where residuals from the ghost points are used
301f746d493SDmitry Karpeev $                 and computed values in ghost regions are added together. Classical
302f746d493SDmitry Karpeev $                 standard additive Schwarz
303f746d493SDmitry Karpeev $  PC_GASM_RESTRICT - residuals from ghost points are used but computed values in ghost
304f746d493SDmitry Karpeev $                    region are discarded. Default
305f746d493SDmitry Karpeev $  PC_GASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
306f746d493SDmitry Karpeev $                       region are added back in
307f746d493SDmitry Karpeev $  PC_GASM_NONE - ghost point residuals are not used, computed ghost values are discarded
308f746d493SDmitry Karpeev $                not very good.
309f746d493SDmitry Karpeev 
310f746d493SDmitry Karpeev    Level: beginner
311f746d493SDmitry Karpeev 
312f746d493SDmitry Karpeev .seealso: PCGASMSetType()
313f746d493SDmitry Karpeev E*/
314f746d493SDmitry Karpeev typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
315f746d493SDmitry Karpeev extern const char *PCGASMTypes[];
316f746d493SDmitry Karpeev 
3177087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
3187087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetTotalSubdomains(PC,PetscInt);
3197087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetOverlap(PC,PetscInt);
3207087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetSortIndices(PC,PetscBool );
321f746d493SDmitry Karpeev 
3227087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetType(PC,PCGASMType);
3237087cfbeSBarry Smith extern PetscErrorCode  PCGASMCreateSubdomains(Mat,PetscInt,IS*[]);
3247087cfbeSBarry Smith extern PetscErrorCode  PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
3257087cfbeSBarry Smith extern PetscErrorCode  PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
3267087cfbeSBarry Smith extern PetscErrorCode  PCGASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
3277087cfbeSBarry Smith extern PetscErrorCode  PCGASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
328f746d493SDmitry Karpeev 
329f746d493SDmitry Karpeev /*E
3303d957683SBarry Smith     PCCompositeType - Determines how two or more preconditioner are composed
3313d957683SBarry Smith 
3323d957683SBarry Smith $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
3333d957683SBarry Smith $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
3343d957683SBarry Smith $                                computed after the previous preconditioner application
335ccb205f8SBarry Smith $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
336ccb205f8SBarry Smith $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
3373d957683SBarry Smith $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
3383d957683SBarry Smith $                         where first preconditioner is built from alpha I + S and second from
3393d957683SBarry Smith $                         alpha I + R
3403d957683SBarry Smith 
3413d957683SBarry Smith    Level: beginner
3423d957683SBarry Smith 
3433d957683SBarry Smith .seealso: PCCompositeSetType()
3443d957683SBarry Smith E*/
3453b224e63SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
3469dcbbd2bSBarry Smith extern const char *PCCompositeTypes[];
3479dcbbd2bSBarry Smith 
3487087cfbeSBarry Smith extern PetscErrorCode  PCCompositeSetUseTrue(PC);
3497087cfbeSBarry Smith extern PetscErrorCode  PCCompositeSetType(PC,PCCompositeType);
3507087cfbeSBarry Smith extern PetscErrorCode  PCCompositeAddPC(PC,PCType);
3517087cfbeSBarry Smith extern PetscErrorCode  PCCompositeGetPC(PC,PetscInt,PC *);
3527087cfbeSBarry Smith extern PetscErrorCode  PCCompositeSpecialSetAlpha(PC,PetscScalar);
353981c4779SBarry Smith 
3547087cfbeSBarry Smith extern PetscErrorCode  PCRedundantSetNumber(PC,PetscInt);
3557087cfbeSBarry Smith extern PetscErrorCode  PCRedundantSetScatter(PC,VecScatter,VecScatter);
3567087cfbeSBarry Smith extern PetscErrorCode  PCRedundantGetOperators(PC,Mat*,Mat*);
357da3a660dSBarry Smith 
3587087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetEpsilon(PC,double);
3597087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetNBSteps(PC,PetscInt);
3607087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetMax(PC,PetscInt);
3617087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetMaxNew(PC,PetscInt);
3627087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetBlockSize(PC,PetscInt);
3637087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetCacheSize(PC,PetscInt);
3647087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetVerbose(PC,PetscInt);
3657087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetSp(PC,PetscInt);
3663304466cSBarry Smith 
3677087cfbeSBarry Smith extern PetscErrorCode  PCHYPRESetType(PC,const char[]);
3687087cfbeSBarry Smith extern PetscErrorCode  PCHYPREGetType(PC,const char*[]);
3697087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
3707087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
3713304466cSBarry Smith 
3727087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*);
3737087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSetType(PC,PCCompositeType);
3747087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSetBlockSize(PC,PetscInt);
3757087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSetIS(PC,const char[],IS);
376084e4875SJed Brown 
377084e4875SJed Brown /*E
378084e4875SJed Brown     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
379084e4875SJed Brown 
380084e4875SJed Brown     Level: intermediate
381084e4875SJed Brown 
382084e4875SJed Brown .seealso: PCFieldSplitSchurPrecondition()
383084e4875SJed Brown E*/
384084e4875SJed Brown typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_DIAG,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
385f7c28744SJed Brown extern const char *const PCFieldSplitSchurPreTypes[];
386084e4875SJed Brown 
3877087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
3887087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
3893d30b1ffSBarry Smith 
3907087cfbeSBarry Smith extern PetscErrorCode  PCGalerkinSetRestriction(PC,Mat);
3917087cfbeSBarry Smith extern PetscErrorCode  PCGalerkinSetInterpolation(PC,Mat);
3922a6744ebSBarry Smith 
3937087cfbeSBarry Smith extern PetscErrorCode  PCSetCoordinates(PC,PetscInt,PetscReal*);
3947087cfbeSBarry Smith extern PetscErrorCode  PCSASetVectors(PC,PetscInt,PetscReal *);
3952102ba4dSSatish Balay 
3967087cfbeSBarry Smith extern PetscErrorCode  PCPythonSetType(PC,const char[]);
39767fac13cSBarry Smith 
3987087cfbeSBarry Smith extern PetscErrorCode  PCSetDM(PC,DM);
3997087cfbeSBarry Smith extern PetscErrorCode  PCGetDM(PC,DM*);
4006c699258SBarry Smith 
4018154be41SBarry Smith extern PetscErrorCode  PCBiCGStabCUSPSetTolerance(PC,PetscReal);
4023d740270SVictor Minden extern PetscErrorCode  PCBiCGStabCUSPSetIterations(PC,PetscInt);
4033d740270SVictor Minden extern PetscErrorCode  PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
4046c699258SBarry Smith 
405e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
406d3b1e0e7SMatthew Knepley 
407123ea438SMatthew Knepley #endif /* __PETSCPC_H */
408