xref: /petsc/include/petscpc.h (revision 7087cfbefd1a42b179f217f9994fb6cb0d0c1824)
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 
10*7087cfbeSBarry 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"
78ee59dc9aSVictor Minden #define PCSACUDA          "sacuda"
7927c67122SBarry Smith #define PCSVD             "svd"
80123ea438SMatthew Knepley 
81123ea438SMatthew Knepley /* Logging support */
82*7087cfbeSBarry Smith extern PetscClassId  PC_CLASSID;
83123ea438SMatthew Knepley 
843d957683SBarry Smith /*E
853d957683SBarry Smith     PCSide - If the preconditioner is to be applied to the left, right
863d957683SBarry Smith      or symmetrically around the operator.
87d03aef70SBarry Smith 
883d957683SBarry Smith    Level: beginner
8921e95762SBarry Smith 
903d957683SBarry Smith .seealso:
913d957683SBarry Smith E*/
92aad2872bSLois Curfman McInnes typedef enum { PC_LEFT,PC_RIGHT,PC_SYMMETRIC } PCSide;
939dcbbd2bSBarry Smith extern const char *PCSides[];
9472b7852fSLois Curfman McInnes 
95*7087cfbeSBarry Smith extern PetscErrorCode  PCCreate(MPI_Comm,PC*);
96*7087cfbeSBarry Smith extern PetscErrorCode  PCSetType(PC,const PCType);
97*7087cfbeSBarry Smith extern PetscErrorCode  PCSetUp(PC);
98*7087cfbeSBarry Smith extern PetscErrorCode  PCSetUpOnBlocks(PC);
99*7087cfbeSBarry Smith extern PetscErrorCode  PCApply(PC,Vec,Vec);
100*7087cfbeSBarry Smith extern PetscErrorCode  PCApplySymmetricLeft(PC,Vec,Vec);
101*7087cfbeSBarry Smith extern PetscErrorCode  PCApplySymmetricRight(PC,Vec,Vec);
102*7087cfbeSBarry Smith extern PetscErrorCode  PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
103*7087cfbeSBarry Smith extern PetscErrorCode  PCApplyTranspose(PC,Vec,Vec);
104*7087cfbeSBarry Smith extern PetscErrorCode  PCApplyTransposeExists(PC,PetscBool *);
105*7087cfbeSBarry Smith extern PetscErrorCode  PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
1064d0a8057SBarry Smith 
1074d0a8057SBarry Smith /*E
1084d0a8057SBarry Smith     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
1094d0a8057SBarry Smith 
1104d0a8057SBarry Smith    Level: advanced
1114d0a8057SBarry Smith 
1124d0a8057SBarry Smith    Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
1134d0a8057SBarry Smith 
1144d0a8057SBarry Smith .seealso: PCApplyRichardson()
1154d0a8057SBarry Smith E*/
1164d0a8057SBarry Smith typedef enum {
1174d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_RTOL               =  2,
1184d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_ATOL               =  3,
1194d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_ITS                =  4,
1204d0a8057SBarry Smith               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
1214d0a8057SBarry Smith 
122*7087cfbeSBarry Smith extern PetscErrorCode  PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
123*7087cfbeSBarry Smith extern PetscErrorCode  PCApplyRichardsonExists(PC,PetscBool *);
124*7087cfbeSBarry Smith extern PetscErrorCode  PCSetInitialGuessNonzero(PC,PetscBool );
12584cb2905SBarry Smith 
126*7087cfbeSBarry Smith extern PetscErrorCode  PCRegisterDestroy(void);
127*7087cfbeSBarry Smith extern PetscErrorCode  PCRegisterAll(const char[]);
128ace3abfcSBarry Smith extern PetscBool  PCRegisterAllCalled;
12984cb2905SBarry Smith 
130*7087cfbeSBarry Smith extern PetscErrorCode  PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
13130de9b25SBarry Smith 
13230de9b25SBarry Smith /*MC
13330de9b25SBarry Smith    PCRegisterDynamic - Adds a method to the preconditioner package.
13430de9b25SBarry Smith 
13530de9b25SBarry Smith    Synopsis:
1361890ba74SBarry Smith    PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC))
13730de9b25SBarry Smith 
13830de9b25SBarry Smith    Not collective
13930de9b25SBarry Smith 
14030de9b25SBarry Smith    Input Parameters:
14130de9b25SBarry Smith +  name_solver - name of a new user-defined solver
14230de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
14330de9b25SBarry Smith .  name_create - name of routine to create method context
14430de9b25SBarry Smith -  routine_create - routine to create method context
14530de9b25SBarry Smith 
14630de9b25SBarry Smith    Notes:
14730de9b25SBarry Smith    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
14830de9b25SBarry Smith 
14930de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
15030de9b25SBarry Smith    is ignored.
15130de9b25SBarry Smith 
15230de9b25SBarry Smith    Sample usage:
15330de9b25SBarry Smith .vb
15430de9b25SBarry Smith    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
15530de9b25SBarry Smith               "MySolverCreate",MySolverCreate);
15630de9b25SBarry Smith .ve
15730de9b25SBarry Smith 
15830de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
15930de9b25SBarry Smith $     PCSetType(pc,"my_solver")
16030de9b25SBarry Smith    or at runtime via the option
16130de9b25SBarry Smith $     -pc_type my_solver
16230de9b25SBarry Smith 
16330de9b25SBarry Smith    Level: advanced
16430de9b25SBarry Smith 
165ab901514SBarry Smith    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
16630de9b25SBarry Smith            occuring in pathname will be replaced with appropriate values.
16730de9b25SBarry Smith          If your function is not being put into a shared library then use PCRegister() instead
16830de9b25SBarry Smith 
16930de9b25SBarry Smith .keywords: PC, register
17030de9b25SBarry Smith 
17130de9b25SBarry Smith .seealso: PCRegisterAll(), PCRegisterDestroy()
17230de9b25SBarry Smith M*/
173aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
174f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
1756df38c32SLois Curfman McInnes #else
176f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
1776df38c32SLois Curfman McInnes #endif
1786df38c32SLois Curfman McInnes 
179*7087cfbeSBarry Smith extern PetscErrorCode  PCDestroy(PC);
180*7087cfbeSBarry Smith extern PetscErrorCode  PCSetFromOptions(PC);
181*7087cfbeSBarry Smith extern PetscErrorCode  PCGetType(PC,const PCType*);
18214c91fddSBarry Smith 
183*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorGetMatrix(PC,Mat*);
184*7087cfbeSBarry Smith extern PetscErrorCode  PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
185*7087cfbeSBarry Smith extern PetscErrorCode  PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
1865b116368SBarry Smith 
187*7087cfbeSBarry Smith extern PetscErrorCode  PCSetOperators(PC,Mat,Mat,MatStructure);
188*7087cfbeSBarry Smith extern PetscErrorCode  PCGetOperators(PC,Mat*,Mat*,MatStructure*);
189*7087cfbeSBarry Smith extern PetscErrorCode  PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
1904b0e389bSBarry Smith 
191*7087cfbeSBarry Smith extern PetscErrorCode  PCView(PC,PetscViewer);
1927bc3d0afSSatish Balay 
193*7087cfbeSBarry Smith extern PetscErrorCode  PCSetOptionsPrefix(PC,const char[]);
194*7087cfbeSBarry Smith extern PetscErrorCode  PCAppendOptionsPrefix(PC,const char[]);
195*7087cfbeSBarry Smith extern PetscErrorCode  PCGetOptionsPrefix(PC,const char*[]);
1968ed539a5SBarry Smith 
197*7087cfbeSBarry Smith extern PetscErrorCode  PCComputeExplicitOperator(PC,Mat*);
19871601f6fSBarry Smith 
199d6913704SBarry Smith /*
200d6913704SBarry Smith       These are used to provide extra scaling of preconditioned
2010f3b3ca1SHong Zhang    operator for time-stepping schemes like in SUNDIALS
202d6913704SBarry Smith */
203*7087cfbeSBarry Smith extern PetscErrorCode  PCGetDiagonalScale(PC,PetscBool *);
204*7087cfbeSBarry Smith extern PetscErrorCode  PCDiagonalScaleLeft(PC,Vec,Vec);
205*7087cfbeSBarry Smith extern PetscErrorCode  PCDiagonalScaleRight(PC,Vec,Vec);
206*7087cfbeSBarry Smith extern PetscErrorCode  PCSetDiagonalScale(PC,Vec);
207d6913704SBarry Smith 
20884cb2905SBarry Smith /* ------------- options specific to particular preconditioners --------- */
209329f5518SBarry Smith 
210*7087cfbeSBarry Smith extern PetscErrorCode  PCJacobiSetUseRowMax(PC);
211*7087cfbeSBarry Smith extern PetscErrorCode  PCJacobiSetUseRowSum(PC);
212*7087cfbeSBarry Smith extern PetscErrorCode  PCJacobiSetUseAbs(PC);
213*7087cfbeSBarry Smith extern PetscErrorCode  PCSORSetSymmetric(PC,MatSORType);
214*7087cfbeSBarry Smith extern PetscErrorCode  PCSORSetOmega(PC,PetscReal);
215*7087cfbeSBarry Smith extern PetscErrorCode  PCSORSetIterations(PC,PetscInt,PetscInt);
216d03aef70SBarry Smith 
217*7087cfbeSBarry Smith extern PetscErrorCode  PCEisenstatSetOmega(PC,PetscReal);
218*7087cfbeSBarry Smith extern PetscErrorCode  PCEisenstatNoDiagonalScaling(PC);
219421c37bdSBarry Smith 
22015aa81f8SBarry Smith #define USE_PRECONDITIONER_MATRIX 0
22115aa81f8SBarry Smith #define USE_TRUE_MATRIX           1
222*7087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiSetUseTrueLocal(PC);
223*7087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
224*7087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
2251eb62cbbSBarry Smith 
226*7087cfbeSBarry Smith extern PetscErrorCode  PCKSPSetUseTrue(PC);
227981c4779SBarry Smith 
228*7087cfbeSBarry Smith extern PetscErrorCode  PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
229*7087cfbeSBarry Smith extern PetscErrorCode  PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
230*7087cfbeSBarry Smith extern PetscErrorCode  PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
231*7087cfbeSBarry Smith extern PetscErrorCode  PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
232*7087cfbeSBarry Smith extern PetscErrorCode  PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
233*7087cfbeSBarry Smith extern PetscErrorCode  PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
234*7087cfbeSBarry Smith extern PetscErrorCode  PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
235*7087cfbeSBarry Smith extern PetscErrorCode  PCShellGetContext(PC,void**);
236*7087cfbeSBarry Smith extern PetscErrorCode  PCShellSetContext(PC,void*);
237*7087cfbeSBarry Smith extern PetscErrorCode  PCShellSetName(PC,const char[]);
238*7087cfbeSBarry Smith extern PetscErrorCode  PCShellGetName(PC,char*[]);
239aabeff55SBarry Smith 
240*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetZeroPivot(PC,PetscReal);
241d90ac83dSHong Zhang 
242*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetShiftType(PC,MatFactorShiftType);
243*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetShiftAmount(PC,PetscReal);
244d90ac83dSHong Zhang 
245*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
246*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
2472401956bSBarry Smith 
248*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetFill(PC,PetscReal);
249*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetColumnPivot(PC,PetscReal);
250*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
251421c37bdSBarry Smith 
252*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetMatOrderingType(PC,const MatOrderingType);
253*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetReuseOrdering(PC,PetscBool );
254*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetReuseFill(PC,PetscBool );
255*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetUseInPlace(PC);
256*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetAllowDiagonalFill(PC);
257*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetPivotInBlocks(PC,PetscBool );
258f5a88c2aSLois Curfman McInnes 
259*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetLevels(PC,PetscInt);
260*7087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
261b35a507dSBarry Smith 
262*7087cfbeSBarry Smith extern PetscErrorCode  PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
263*7087cfbeSBarry Smith extern PetscErrorCode  PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
264*7087cfbeSBarry Smith extern PetscErrorCode  PCASMSetOverlap(PC,PetscInt);
265*7087cfbeSBarry Smith extern PetscErrorCode  PCASMSetSortIndices(PC,PetscBool );
266f746d493SDmitry Karpeev 
2673d957683SBarry Smith /*E
2683d957683SBarry Smith     PCASMType - Type of additive Schwarz method to use
2693d957683SBarry Smith 
2703d957683SBarry Smith $  PC_ASM_BASIC - symmetric version where residuals from the ghost points are used
2713d957683SBarry Smith $                 and computed values in ghost regions are added together. Classical
2723d957683SBarry Smith $                 standard additive Schwarz
2733d957683SBarry Smith $  PC_ASM_RESTRICT - residuals from ghost points are used but computed values in ghost
2743d957683SBarry Smith $                    region are discarded. Default
2753d957683SBarry Smith $  PC_ASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
2763d957683SBarry Smith $                       region are added back in
2773d957683SBarry Smith $  PC_ASM_NONE - ghost point residuals are not used, computed ghost values are discarded
2783d957683SBarry Smith $                not very good.
2793d957683SBarry Smith 
2803d957683SBarry Smith    Level: beginner
2813d957683SBarry Smith 
2823d957683SBarry Smith .seealso: PCASMSetType()
2833d957683SBarry Smith E*/
284d252947aSBarry Smith typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
2859dcbbd2bSBarry Smith extern const char *PCASMTypes[];
2869dcbbd2bSBarry Smith 
287*7087cfbeSBarry Smith extern PetscErrorCode  PCASMSetType(PC,PCASMType);
288*7087cfbeSBarry Smith extern PetscErrorCode  PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
289*7087cfbeSBarry Smith extern PetscErrorCode  PCASMDestroySubdomains(PetscInt,IS[],IS[]);
290*7087cfbeSBarry Smith extern PetscErrorCode  PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
291*7087cfbeSBarry Smith extern PetscErrorCode  PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
292*7087cfbeSBarry Smith extern PetscErrorCode  PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
293981c4779SBarry Smith 
2943d957683SBarry Smith /*E
295f746d493SDmitry Karpeev     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per domain)
296f746d493SDmitry Karpeev 
297f746d493SDmitry Karpeev $  PC_GASM_BASIC - symmetric version where residuals from the ghost points are used
298f746d493SDmitry Karpeev $                 and computed values in ghost regions are added together. Classical
299f746d493SDmitry Karpeev $                 standard additive Schwarz
300f746d493SDmitry Karpeev $  PC_GASM_RESTRICT - residuals from ghost points are used but computed values in ghost
301f746d493SDmitry Karpeev $                    region are discarded. Default
302f746d493SDmitry Karpeev $  PC_GASM_INTERPOLATE - residuals from ghost points are not used, computed values in ghost
303f746d493SDmitry Karpeev $                       region are added back in
304f746d493SDmitry Karpeev $  PC_GASM_NONE - ghost point residuals are not used, computed ghost values are discarded
305f746d493SDmitry Karpeev $                not very good.
306f746d493SDmitry Karpeev 
307f746d493SDmitry Karpeev    Level: beginner
308f746d493SDmitry Karpeev 
309f746d493SDmitry Karpeev .seealso: PCGASMSetType()
310f746d493SDmitry Karpeev E*/
311f746d493SDmitry Karpeev typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
312f746d493SDmitry Karpeev extern const char *PCGASMTypes[];
313f746d493SDmitry Karpeev 
314*7087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
315*7087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetTotalSubdomains(PC,PetscInt);
316*7087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetOverlap(PC,PetscInt);
317*7087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetSortIndices(PC,PetscBool );
318f746d493SDmitry Karpeev 
319*7087cfbeSBarry Smith extern PetscErrorCode  PCGASMSetType(PC,PCGASMType);
320*7087cfbeSBarry Smith extern PetscErrorCode  PCGASMCreateSubdomains(Mat,PetscInt,IS*[]);
321*7087cfbeSBarry Smith extern PetscErrorCode  PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
322*7087cfbeSBarry Smith extern PetscErrorCode  PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
323*7087cfbeSBarry Smith extern PetscErrorCode  PCGASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
324*7087cfbeSBarry Smith extern PetscErrorCode  PCGASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
325f746d493SDmitry Karpeev 
326f746d493SDmitry Karpeev /*E
3273d957683SBarry Smith     PCCompositeType - Determines how two or more preconditioner are composed
3283d957683SBarry Smith 
3293d957683SBarry Smith $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
3303d957683SBarry Smith $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
3313d957683SBarry Smith $                                computed after the previous preconditioner application
332ccb205f8SBarry Smith $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
333ccb205f8SBarry Smith $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
3343d957683SBarry Smith $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
3353d957683SBarry Smith $                         where first preconditioner is built from alpha I + S and second from
3363d957683SBarry Smith $                         alpha I + R
3373d957683SBarry Smith 
3383d957683SBarry Smith    Level: beginner
3393d957683SBarry Smith 
3403d957683SBarry Smith .seealso: PCCompositeSetType()
3413d957683SBarry Smith E*/
3423b224e63SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
3439dcbbd2bSBarry Smith extern const char *PCCompositeTypes[];
3449dcbbd2bSBarry Smith 
345*7087cfbeSBarry Smith extern PetscErrorCode  PCCompositeSetUseTrue(PC);
346*7087cfbeSBarry Smith extern PetscErrorCode  PCCompositeSetType(PC,PCCompositeType);
347*7087cfbeSBarry Smith extern PetscErrorCode  PCCompositeAddPC(PC,PCType);
348*7087cfbeSBarry Smith extern PetscErrorCode  PCCompositeGetPC(PC,PetscInt,PC *);
349*7087cfbeSBarry Smith extern PetscErrorCode  PCCompositeSpecialSetAlpha(PC,PetscScalar);
350981c4779SBarry Smith 
351*7087cfbeSBarry Smith extern PetscErrorCode  PCRedundantSetNumber(PC,PetscInt);
352*7087cfbeSBarry Smith extern PetscErrorCode  PCRedundantSetScatter(PC,VecScatter,VecScatter);
353*7087cfbeSBarry Smith extern PetscErrorCode  PCRedundantGetOperators(PC,Mat*,Mat*);
354*7087cfbeSBarry Smith extern PetscErrorCode  PCRedundantGetPC(PC,PC*);
355da3a660dSBarry Smith 
356*7087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetEpsilon(PC,double);
357*7087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetNBSteps(PC,PetscInt);
358*7087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetMax(PC,PetscInt);
359*7087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetMaxNew(PC,PetscInt);
360*7087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetBlockSize(PC,PetscInt);
361*7087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetCacheSize(PC,PetscInt);
362*7087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetVerbose(PC,PetscInt);
363*7087cfbeSBarry Smith extern PetscErrorCode  PCSPAISetSp(PC,PetscInt);
3643304466cSBarry Smith 
365*7087cfbeSBarry Smith extern PetscErrorCode  PCHYPRESetType(PC,const char[]);
366*7087cfbeSBarry Smith extern PetscErrorCode  PCHYPREGetType(PC,const char*[]);
367*7087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
368*7087cfbeSBarry Smith extern PetscErrorCode  PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
3693304466cSBarry Smith 
370*7087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*);
371*7087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSetType(PC,PCCompositeType);
372*7087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSetBlockSize(PC,PetscInt);
373*7087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSetIS(PC,const char[],IS);
374084e4875SJed Brown 
375084e4875SJed Brown /*E
376084e4875SJed Brown     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
377084e4875SJed Brown 
378084e4875SJed Brown     Level: intermediate
379084e4875SJed Brown 
380084e4875SJed Brown .seealso: PCFieldSplitSchurPrecondition()
381084e4875SJed Brown E*/
382084e4875SJed Brown typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_DIAG,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
383f7c28744SJed Brown extern const char *const PCFieldSplitSchurPreTypes[];
384084e4875SJed Brown 
385*7087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
386*7087cfbeSBarry Smith extern PetscErrorCode  PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
3873d30b1ffSBarry Smith 
388*7087cfbeSBarry Smith extern PetscErrorCode  PCGalerkinSetRestriction(PC,Mat);
389*7087cfbeSBarry Smith extern PetscErrorCode  PCGalerkinSetInterpolation(PC,Mat);
3902a6744ebSBarry Smith 
391*7087cfbeSBarry Smith extern PetscErrorCode  PCSetCoordinates(PC,PetscInt,PetscReal*);
392*7087cfbeSBarry Smith extern PetscErrorCode  PCSASetVectors(PC,PetscInt,PetscReal *);
3932102ba4dSSatish Balay 
394*7087cfbeSBarry Smith extern PetscErrorCode  PCPythonSetType(PC,const char[]);
39567fac13cSBarry Smith 
396*7087cfbeSBarry Smith extern PetscErrorCode  PCSetDM(PC,DM);
397*7087cfbeSBarry Smith extern PetscErrorCode  PCGetDM(PC,DM*);
3986c699258SBarry Smith 
3996c699258SBarry Smith 
400e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
401d3b1e0e7SMatthew Knepley 
402123ea438SMatthew Knepley #endif /* __PETSCPC_H */
403