xref: /petsc/include/petscpc.h (revision 1a83f524638a4da8317a6bd80eb6d9a2936d8384)
1d03aef70SBarry Smith /*
237f753daSBarry Smith       Preconditioner module.
3d03aef70SBarry Smith */
40a835dfdSSatish Balay #if !defined(__PETSCPC_H)
50a835dfdSSatish Balay #define __PETSCPC_H
62c8e378dSBarry Smith #include <petscdm.h>
7d03aef70SBarry Smith 
8014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCInitializePackage(const char[]);
91dbb0a54SBarry Smith 
10eec0b4cfSBarry Smith /*
11eec0b4cfSBarry Smith     PCList contains the list of preconditioners currently registered
12f1af5d2fSBarry Smith    These are added with the PCRegisterDynamic() macro
13eec0b4cfSBarry Smith */
14014dd563SJed Brown PETSC_EXTERN PetscFList PCList;
1582bf6240SBarry Smith 
163d957683SBarry Smith /*S
173d957683SBarry Smith      PC - Abstract PETSc object that manages all preconditioners
183d957683SBarry Smith 
193d957683SBarry Smith    Level: beginner
203d957683SBarry Smith 
213d957683SBarry Smith   Concepts: preconditioners
223d957683SBarry Smith 
231a480d89SAdministrator .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
243d957683SBarry Smith S*/
253d957683SBarry Smith typedef struct _p_PC* PC;
263d957683SBarry Smith 
2776bdecfbSBarry Smith /*J
283d957683SBarry Smith     PCType - String with the name of a PETSc preconditioner method or the creation function
293d957683SBarry Smith        with an optional dynamic library name, for example
303d957683SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mypccreate()
313d957683SBarry Smith 
323d957683SBarry Smith    Level: beginner
333d957683SBarry Smith 
341a480d89SAdministrator    Notes: Click on the links below to see details on a particular solver
351a480d89SAdministrator 
361a480d89SAdministrator .seealso: PCSetType(), PC, PCCreate()
3776bdecfbSBarry Smith J*/
3819fd82e9SBarry Smith typedef const char* PCType;
3982bf6240SBarry Smith #define PCNONE            "none"
4082bf6240SBarry Smith #define PCJACOBI          "jacobi"
4182bf6240SBarry Smith #define PCSOR             "sor"
4282bf6240SBarry Smith #define PCLU              "lu"
4382bf6240SBarry Smith #define PCSHELL           "shell"
4482bf6240SBarry Smith #define PCBJACOBI         "bjacobi"
4582bf6240SBarry Smith #define PCMG              "mg"
4682bf6240SBarry Smith #define PCEISENSTAT       "eisenstat"
4782bf6240SBarry Smith #define PCILU             "ilu"
4882bf6240SBarry Smith #define PCICC             "icc"
4982bf6240SBarry Smith #define PCASM             "asm"
50ab3e923fSDmitry Karpeev #define PCGASM            "gasm"
5194b7f48cSBarry Smith #define PCKSP             "ksp"
5282bf6240SBarry Smith #define PCCOMPOSITE       "composite"
53421c37bdSBarry Smith #define PCREDUNDANT       "redundant"
5427b520f0SBarry Smith #define PCSPAI            "spai"
55186905e3SBarry Smith #define PCNN              "nn"
564bbc92c1SBarry Smith #define PCCHOLESKY        "cholesky"
573a7fca6bSBarry Smith #define PCPBJACOBI        "pbjacobi"
587f5ff6fdSBarry Smith #define PCMAT             "mat"
59c4888f26SBarry Smith #define PCHYPRE           "hypre"
6037f80224SJose E. Roman #define PCPARMS           "parms"
610971522cSBarry Smith #define PCFIELDSPLIT      "fieldsplit"
62be16f70fSSatish Balay #define PCTFS             "tfs"
635582bec1SHong Zhang #define PCML              "ml"
642a6744ebSBarry Smith #define PCGALERKIN        "galerkin"
657233f9f0SBarry Smith #define PCEXOTIC          "exotic"
66bad7cb1dSBarry Smith #define PCHMPI            "hmpi"
67cf037197Ssdaitch #define PCSUPPORTGRAPH    "supportgraph"
68f4b8409dSBarry Smith #define PCASA             "asa"
6924c02a0fSBarry Smith #define PCCP              "cp"
70628b8437SMatthew Knepley #define PCBFBT            "bfbt"
71519d70e2SJed Brown #define PCLSC             "lsc"
721d6018f0SLisandro Dalcin #define PCPYTHON          "python"
73f91d8e95SBarry Smith #define PCPFMG            "pfmg"
74d851a50bSGlenn Hammond #define PCSYSPFMG         "syspfmg"
75df826632SBarry Smith #define PCREDISTRIBUTE    "redistribute"
76ac793be5SBarry Smith #define PCSVD             "svd"
77ac793be5SBarry Smith #define PCGAMG            "gamg"
78ac793be5SBarry Smith #define PCSACUSP          "sacusp"        /* these four run on NVIDIA GPUs using CUSP */
79bfc29b71SVictor Minden #define PCSACUSPPOLY      "sacusppoly"
808154be41SBarry Smith #define PCBICGSTABCUSP    "bicgstabcusp"
8104b59e88SVictor Minden #define PCAINVCUSP        "ainvcusp"
820c7d97c5SJed Brown #define PCBDDC            "bddc"
83123ea438SMatthew Knepley 
84123ea438SMatthew Knepley /* Logging support */
85014dd563SJed Brown PETSC_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*/
959e568555SJed Brown typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
969e568555SJed Brown #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
970910f4e4SJed Brown PETSC_EXTERN const char *const *const PCSides;
9872b7852fSLois Curfman McInnes 
99014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
10019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType);
101014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetUp(PC);
102014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
103014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
104014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
106014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *);
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
1104d0a8057SBarry Smith 
1114d0a8057SBarry Smith /*E
1124d0a8057SBarry Smith     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
1134d0a8057SBarry Smith 
1144d0a8057SBarry Smith    Level: advanced
1154d0a8057SBarry Smith 
1164d0a8057SBarry Smith    Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
1174d0a8057SBarry Smith 
1184d0a8057SBarry Smith .seealso: PCApplyRichardson()
1194d0a8057SBarry Smith E*/
1204d0a8057SBarry Smith typedef enum {
1214d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_RTOL               =  2,
1224d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_ATOL               =  3,
1234d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_ITS                =  4,
1244d0a8057SBarry Smith               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
1254d0a8057SBarry Smith 
126014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *);
128014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool );
12984cb2905SBarry Smith 
130014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRegisterDestroy(void);
131014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRegisterAll(const char[]);
132014dd563SJed Brown PETSC_EXTERN PetscBool PCRegisterAllCalled;
13384cb2905SBarry Smith 
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
13530de9b25SBarry Smith 
13630de9b25SBarry Smith /*MC
13730de9b25SBarry Smith    PCRegisterDynamic - Adds a method to the preconditioner package.
13830de9b25SBarry Smith 
13930de9b25SBarry Smith    Synopsis:
1401890ba74SBarry Smith    PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC))
14130de9b25SBarry Smith 
14230de9b25SBarry Smith    Not collective
14330de9b25SBarry Smith 
14430de9b25SBarry Smith    Input Parameters:
14530de9b25SBarry Smith +  name_solver - name of a new user-defined solver
14630de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
14730de9b25SBarry Smith .  name_create - name of routine to create method context
14830de9b25SBarry Smith -  routine_create - routine to create method context
14930de9b25SBarry Smith 
15030de9b25SBarry Smith    Notes:
15130de9b25SBarry Smith    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
15230de9b25SBarry Smith 
15330de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
15430de9b25SBarry Smith    is ignored.
15530de9b25SBarry Smith 
15630de9b25SBarry Smith    Sample usage:
15730de9b25SBarry Smith .vb
15830de9b25SBarry Smith    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
15930de9b25SBarry Smith               "MySolverCreate",MySolverCreate);
16030de9b25SBarry Smith .ve
16130de9b25SBarry Smith 
16230de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
16330de9b25SBarry Smith $     PCSetType(pc,"my_solver")
16430de9b25SBarry Smith    or at runtime via the option
16530de9b25SBarry Smith $     -pc_type my_solver
16630de9b25SBarry Smith 
16730de9b25SBarry Smith    Level: advanced
16830de9b25SBarry Smith 
169ab901514SBarry Smith    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
17030de9b25SBarry Smith            occuring in pathname will be replaced with appropriate values.
17130de9b25SBarry Smith          If your function is not being put into a shared library then use PCRegister() instead
17230de9b25SBarry Smith 
17330de9b25SBarry Smith .keywords: PC, register
17430de9b25SBarry Smith 
17530de9b25SBarry Smith .seealso: PCRegisterAll(), PCRegisterDestroy()
17630de9b25SBarry Smith M*/
177aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
178f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
1796df38c32SLois Curfman McInnes #else
180f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
1816df38c32SLois Curfman McInnes #endif
1826df38c32SLois Curfman McInnes 
183014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCReset(PC);
184014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
185014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);
18619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*);
18714c91fddSBarry Smith 
188014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
189014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
190014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
1915b116368SBarry Smith 
192014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
193014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
194014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
1954b0e389bSBarry Smith 
196014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer);
1977bc3d0afSSatish Balay 
198014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
199014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
200014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);
2018ed539a5SBarry Smith 
202014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
20371601f6fSBarry Smith 
204d6913704SBarry Smith /*
205d6913704SBarry Smith       These are used to provide extra scaling of preconditioned
2060f3b3ca1SHong Zhang    operator for time-stepping schemes like in SUNDIALS
207d6913704SBarry Smith */
208014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *);
209014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
210014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
211014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec);
212d6913704SBarry Smith 
21384cb2905SBarry Smith /* ------------- options specific to particular preconditioners --------- */
214329f5518SBarry Smith 
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC);
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC);
218014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
219014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
220014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
221d03aef70SBarry Smith 
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
224421c37bdSBarry Smith 
22515aa81f8SBarry Smith #define USE_PRECONDITIONER_MATRIX 0
22615aa81f8SBarry Smith #define USE_TRUE_MATRIX           1
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC);
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
2301eb62cbbSBarry Smith 
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPSetUseTrue(PC);
232981c4779SBarry Smith 
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
234014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
237014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
238014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
239014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
240014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**);
241014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*);
242014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
243014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]);
244aabeff55SBarry Smith 
245014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);
246d90ac83dSHong Zhang 
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType);
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal);
249d90ac83dSHong Zhang 
250014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC);
2532401956bSBarry Smith 
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
257421c37bdSBarry Smith 
25819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType);
259014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool );
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool );
261014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC);
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool );
264f5a88c2aSLois Curfman McInnes 
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
267b35a507dSBarry Smith 
268014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
269014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
271014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool );
272f746d493SDmitry Karpeev 
2733d957683SBarry Smith /*E
2743d957683SBarry Smith     PCASMType - Type of additive Schwarz method to use
2753d957683SBarry Smith 
276feb221f8SDmitry Karpeev $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
277feb221f8SDmitry Karpeev $                        and computed values in ghost regions are added together.
278feb221f8SDmitry Karpeev $                        Classical standard additive Schwarz.
279feb221f8SDmitry Karpeev $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
280feb221f8SDmitry Karpeev $                        region are discarded.
281feb221f8SDmitry Karpeev $                        Default.
282feb221f8SDmitry Karpeev $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
283feb221f8SDmitry Karpeev $                        region are added back in.
284feb221f8SDmitry Karpeev $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
285feb221f8SDmitry Karpeev $                        discarded.
286feb221f8SDmitry Karpeev $                        Not very good.
2873d957683SBarry Smith 
2883d957683SBarry Smith    Level: beginner
2893d957683SBarry Smith 
2903d957683SBarry Smith .seealso: PCASMSetType()
2913d957683SBarry Smith E*/
292d252947aSBarry Smith typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
2936a6fc655SJed Brown PETSC_EXTERN const char *const PCASMTypes[];
2949dcbbd2bSBarry Smith 
295014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
296014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
297014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
298014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
299014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
300014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
301981c4779SBarry Smith 
3023d957683SBarry Smith /*E
3036a4f0f73SDmitry Karpeev     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
304f746d493SDmitry Karpeev 
3056a4f0f73SDmitry Karpeev    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
3066a4f0f73SDmitry Karpeev    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
3076a4f0f73SDmitry Karpeev    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
3086a4f0f73SDmitry Karpeev    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
3096a4f0f73SDmitry Karpeev 
3106a4f0f73SDmitry Karpeev $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
3116a4f0f73SDmitry Karpeev $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
3126a4f0f73SDmitry Karpeev $                        from neighboring subdomains.
313feb221f8SDmitry Karpeev $                        Classical standard additive Schwarz.
3146a4f0f73SDmitry Karpeev $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
3156a4f0f73SDmitry Karpeev $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
3166a4f0f73SDmitry Karpeev $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
3176a4f0f73SDmitry Karpeev $                        assumption).
318feb221f8SDmitry Karpeev $                        Default.
3196a4f0f73SDmitry Karpeev $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
3206a4f0f73SDmitry Karpeev $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
3216a4f0f73SDmitry Karpeev $                        from neighboring subdomains.
3226a4f0f73SDmitry Karpeev $
3236a4f0f73SDmitry Karpeev $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
324feb221f8SDmitry Karpeev $                        Not very good.
325f746d493SDmitry Karpeev 
326f746d493SDmitry Karpeev    Level: beginner
327f746d493SDmitry Karpeev 
328f746d493SDmitry Karpeev .seealso: PCGASMSetType()
329f746d493SDmitry Karpeev E*/
330f746d493SDmitry Karpeev typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
3316a6fc655SJed Brown PETSC_EXTERN const char *const PCGASMTypes[];
332f746d493SDmitry Karpeev 
333014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
334014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool);
335014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
336014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );
337f746d493SDmitry Karpeev 
338014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
339014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]);
340014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
341014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
342014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
343014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);
344f746d493SDmitry Karpeev 
345f746d493SDmitry Karpeev /*E
3463d957683SBarry Smith     PCCompositeType - Determines how two or more preconditioner are composed
3473d957683SBarry Smith 
3483d957683SBarry Smith $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
3493d957683SBarry Smith $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
3503d957683SBarry Smith $                                computed after the previous preconditioner application
351ccb205f8SBarry Smith $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
352ccb205f8SBarry Smith $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
3533d957683SBarry Smith $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
3543d957683SBarry Smith $                         where first preconditioner is built from alpha I + S and second from
3553d957683SBarry Smith $                         alpha I + R
3563d957683SBarry Smith 
3573d957683SBarry Smith    Level: beginner
3583d957683SBarry Smith 
3593d957683SBarry Smith .seealso: PCCompositeSetType()
3603d957683SBarry Smith E*/
3613b224e63SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
3626a6fc655SJed Brown PETSC_EXTERN const char *const PCCompositeTypes[];
3639dcbbd2bSBarry Smith 
364014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSetUseTrue(PC);
365014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
36619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
367014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
368014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
369981c4779SBarry Smith 
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
373da3a660dSBarry Smith 
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
376014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
377014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
378014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
379014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
380014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
381014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
3823304466cSBarry Smith 
383014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
384014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
385014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
386014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
3873304466cSBarry Smith 
388014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
389014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
390014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
391014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
392014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
393014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
394084e4875SJed Brown 
395084e4875SJed Brown /*E
396084e4875SJed Brown     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
397084e4875SJed Brown 
398084e4875SJed Brown     Level: intermediate
399084e4875SJed Brown 
400084e4875SJed Brown .seealso: PCFieldSplitSchurPrecondition()
401084e4875SJed Brown E*/
402084e4875SJed Brown typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_DIAG,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
403014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
404084e4875SJed Brown 
405ab1df9f5SJed Brown /*E
406c9c6ffaaSJed Brown     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
407ab1df9f5SJed Brown 
408ab1df9f5SJed Brown     Level: intermediate
409ab1df9f5SJed Brown 
410c9c6ffaaSJed Brown .seealso: PCFieldSplitSetSchurFactType()
411ab1df9f5SJed Brown E*/
412ab1df9f5SJed Brown typedef enum {
413c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
414c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
415c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
416c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_FULL
417c9c6ffaaSJed Brown } PCFieldSplitSchurFactType;
418014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
419ab1df9f5SJed Brown 
420014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
421014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
422014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
4233d30b1ffSBarry Smith 
424014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
425014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
4262a6744ebSBarry Smith 
427014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
428014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);
4292102ba4dSSatish Balay 
430014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);
43167fac13cSBarry Smith 
432014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
433014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);
4346c699258SBarry Smith 
435014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
436014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);
4371b2093e4SBarry Smith 
438014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
439014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
440014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
4416c699258SBarry Smith 
442014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
443014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
444014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
445014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
44637f80224SJose E. Roman /*E
44737f80224SJose E. Roman     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
44837f80224SJose E. Roman 
44937f80224SJose E. Roman     Level: intermediate
45037f80224SJose E. Roman 
45137f80224SJose E. Roman .seealso: PCPARMSSetGlobal()
45237f80224SJose E. Roman E*/
45337f80224SJose E. Roman typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
4546a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
45537f80224SJose E. Roman /*E
45637f80224SJose E. Roman     PCPARMSLocalType - Determines the local preconditioner method in PARMS
45737f80224SJose E. Roman 
45837f80224SJose E. Roman     Level: intermediate
45937f80224SJose E. Roman 
46037f80224SJose E. Roman .seealso: PCPARMSSetLocal()
46137f80224SJose E. Roman E*/
46237f80224SJose E. Roman typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
4636a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSLocalTypes[];
46437f80224SJose E. Roman 
465014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type);
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type);
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits);
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart);
469014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym);
470014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2);
47137f80224SJose E. Roman 
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
473014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
474014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
475014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
476014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
478014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
47919fd82e9SBarry Smith typedef const char* PCGAMGType;
48019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType );
481014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n);
482014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n);
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool);
48456de70efSSatish Balay 
4850c7d97c5SJed Brown #if defined(PETSC_HAVE_PCBDDC)
48653cdbc3dSStefano Zampini /* Enum defining how to treat the coarse problem */
4870c7d97c5SJed Brown typedef enum {SEQUENTIAL_BDDC,REPLICATED_BDDC,PARALLEL_BDDC,MULTILEVEL_BDDC} CoarseProblemType;
4884fad6a16SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
4894fad6a16SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetMaxLevels(PC,PetscInt);
4900bdf917eSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace);
491014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
495014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseProblemType(PC,CoarseProblemType);
496014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
497*1a83f524SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode);
4983425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*);
4993425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec);
5003425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec);
5010c7d97c5SJed Brown #endif
5020c7d97c5SJed Brown 
503b965d45eSStefano Zampini PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool);
504014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
50546caae47SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec);
506d3b1e0e7SMatthew Knepley 
507123ea438SMatthew Knepley #endif /* __PETSCPC_H */
508