xref: /petsc/include/petscpc.h (revision 49517cde1b273295cf3941ef5a8aecdb444b7da7)
1d03aef70SBarry Smith /*
237f753daSBarry Smith       Preconditioner module.
3d03aef70SBarry Smith */
40a835dfdSSatish Balay #if !defined(__PETSCPC_H)
50a835dfdSSatish Balay #define __PETSCPC_H
61e25c274SJed Brown #include <petscmat.h>
71e25c274SJed Brown #include <petscdmtypes.h>
8d03aef70SBarry Smith 
9014dd563SJed Brown PETSC_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 */
15140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList 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 
2876bdecfbSBarry Smith /*J
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()
3876bdecfbSBarry Smith J*/
3919fd82e9SBarry Smith typedef const char* PCType;
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"
6137f80224SJose E. Roman #define PCPARMS           "parms"
620971522cSBarry Smith #define PCFIELDSPLIT      "fieldsplit"
63be16f70fSSatish Balay #define PCTFS             "tfs"
645582bec1SHong Zhang #define PCML              "ml"
652a6744ebSBarry Smith #define PCGALERKIN        "galerkin"
667233f9f0SBarry Smith #define PCEXOTIC          "exotic"
67bad7cb1dSBarry Smith #define PCHMPI            "hmpi"
68cf037197Ssdaitch #define PCSUPPORTGRAPH    "supportgraph"
69f4b8409dSBarry Smith #define PCASA             "asa"
7024c02a0fSBarry Smith #define PCCP              "cp"
71628b8437SMatthew Knepley #define PCBFBT            "bfbt"
72519d70e2SJed Brown #define PCLSC             "lsc"
731d6018f0SLisandro Dalcin #define PCPYTHON          "python"
74f91d8e95SBarry Smith #define PCPFMG            "pfmg"
75d851a50bSGlenn Hammond #define PCSYSPFMG         "syspfmg"
76df826632SBarry Smith #define PCREDISTRIBUTE    "redistribute"
77ac793be5SBarry Smith #define PCSVD             "svd"
78ac793be5SBarry Smith #define PCGAMG            "gamg"
79ac793be5SBarry Smith #define PCSACUSP          "sacusp"        /* these four run on NVIDIA GPUs using CUSP */
80bfc29b71SVictor Minden #define PCSACUSPPOLY      "sacusppoly"
818154be41SBarry Smith #define PCBICGSTABCUSP    "bicgstabcusp"
8204b59e88SVictor Minden #define PCAINVCUSP        "ainvcusp"
830c7d97c5SJed Brown #define PCBDDC            "bddc"
84123ea438SMatthew Knepley 
85123ea438SMatthew Knepley /* Logging support */
86014dd563SJed Brown PETSC_EXTERN PetscClassId PC_CLASSID;
87123ea438SMatthew Knepley 
883d957683SBarry Smith /*E
893d957683SBarry Smith     PCSide - If the preconditioner is to be applied to the left, right
903d957683SBarry Smith      or symmetrically around the operator.
91d03aef70SBarry Smith 
923d957683SBarry Smith    Level: beginner
9321e95762SBarry Smith 
943d957683SBarry Smith .seealso:
953d957683SBarry Smith E*/
969e568555SJed Brown typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
979e568555SJed Brown #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
980910f4e4SJed Brown PETSC_EXTERN const char *const *const PCSides;
9972b7852fSLois Curfman McInnes 
100014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*);
10119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType);
102014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetUp(PC);
103014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC);
104014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec);
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec);
106014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec);
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec);
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec);
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *);
110014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec);
1114d0a8057SBarry Smith 
11255849f57SBarry Smith #define PC_FILE_CLASSID 1211222
11355849f57SBarry Smith 
1144d0a8057SBarry Smith /*E
1154d0a8057SBarry Smith     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
1164d0a8057SBarry Smith 
1174d0a8057SBarry Smith    Level: advanced
1184d0a8057SBarry Smith 
1194d0a8057SBarry Smith    Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
1204d0a8057SBarry Smith 
1214d0a8057SBarry Smith .seealso: PCApplyRichardson()
1224d0a8057SBarry Smith E*/
1234d0a8057SBarry Smith typedef enum {
1244d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_RTOL               =  2,
1254d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_ATOL               =  3,
1264d0a8057SBarry Smith               PCRICHARDSON_CONVERGED_ITS                =  4,
1274d0a8057SBarry Smith               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
1284d0a8057SBarry Smith 
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*);
130014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *);
131014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool );
132*49517cdeSBarry Smith PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool);
133*49517cdeSBarry Smith PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*);
13484cb2905SBarry Smith 
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRegisterDestroy(void);
136014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRegisterAll(const char[]);
137014dd563SJed Brown PETSC_EXTERN PetscBool PCRegisterAllCalled;
13884cb2905SBarry Smith 
139014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC));
14030de9b25SBarry Smith 
14130de9b25SBarry Smith /*MC
14230de9b25SBarry Smith    PCRegisterDynamic - Adds a method to the preconditioner package.
14330de9b25SBarry Smith 
14430de9b25SBarry Smith    Synopsis:
145f2ba6396SBarry Smith     #include "petscpc.h"
1461890ba74SBarry Smith    PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC))
14730de9b25SBarry Smith 
14830de9b25SBarry Smith    Not collective
14930de9b25SBarry Smith 
15030de9b25SBarry Smith    Input Parameters:
15130de9b25SBarry Smith +  name_solver - name of a new user-defined solver
15230de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
15330de9b25SBarry Smith .  name_create - name of routine to create method context
15430de9b25SBarry Smith -  routine_create - routine to create method context
15530de9b25SBarry Smith 
15630de9b25SBarry Smith    Notes:
15730de9b25SBarry Smith    PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners.
15830de9b25SBarry Smith 
15930de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
16030de9b25SBarry Smith    is ignored.
16130de9b25SBarry Smith 
16230de9b25SBarry Smith    Sample usage:
16330de9b25SBarry Smith .vb
16430de9b25SBarry Smith    PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib",
16530de9b25SBarry Smith               "MySolverCreate",MySolverCreate);
16630de9b25SBarry Smith .ve
16730de9b25SBarry Smith 
16830de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
16930de9b25SBarry Smith $     PCSetType(pc,"my_solver")
17030de9b25SBarry Smith    or at runtime via the option
17130de9b25SBarry Smith $     -pc_type my_solver
17230de9b25SBarry Smith 
17330de9b25SBarry Smith    Level: advanced
17430de9b25SBarry Smith 
175ab901514SBarry Smith    Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},  or ${any environmental variable}
17630de9b25SBarry Smith            occuring in pathname will be replaced with appropriate values.
17730de9b25SBarry Smith          If your function is not being put into a shared library then use PCRegister() instead
17830de9b25SBarry Smith 
17930de9b25SBarry Smith .keywords: PC, register
18030de9b25SBarry Smith 
18130de9b25SBarry Smith .seealso: PCRegisterAll(), PCRegisterDestroy()
18230de9b25SBarry Smith M*/
183aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
184f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0)
1856df38c32SLois Curfman McInnes #else
186f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d)
1876df38c32SLois Curfman McInnes #endif
1886df38c32SLois Curfman McInnes 
189014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCReset(PC);
190014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDestroy(PC*);
191014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC);
19219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*);
19314c91fddSBarry Smith 
194014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*);
195014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*);
196014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*);
1975b116368SBarry Smith 
198014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure);
199014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*);
200014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *);
2014b0e389bSBarry Smith 
202014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer);
20355849f57SBarry Smith PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer);
2047bc3d0afSSatish Balay 
205014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]);
206014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]);
207014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]);
2088ed539a5SBarry Smith 
209014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*);
21071601f6fSBarry Smith 
211d6913704SBarry Smith /*
212d6913704SBarry Smith       These are used to provide extra scaling of preconditioned
2130f3b3ca1SHong Zhang    operator for time-stepping schemes like in SUNDIALS
214d6913704SBarry Smith */
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *);
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec);
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec);
218014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec);
219d6913704SBarry Smith 
22084cb2905SBarry Smith /* ------------- options specific to particular preconditioners --------- */
221329f5518SBarry Smith 
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC);
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC);
224014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC);
225014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType);
226014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal);
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt);
228d03aef70SBarry Smith 
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal);
230014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC);
231421c37bdSBarry Smith 
232014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]);
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]);
2341eb62cbbSBarry Smith 
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec));
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec));
237014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec));
238014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC));
239014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*));
240014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer));
241014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC));
242014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**);
243014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*);
244014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]);
245014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]);
246aabeff55SBarry Smith 
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal);
248d90ac83dSHong Zhang 
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType);
250014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal);
251d90ac83dSHong Zhang 
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage);
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*);
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC);
2552401956bSBarry Smith 
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal);
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal);
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal);
259421c37bdSBarry Smith 
26019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType);
261014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool );
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool );
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool );
266f5a88c2aSLois Curfman McInnes 
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt);
268014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt);
269b35a507dSBarry Smith 
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]);
271014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]);
272014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt);
273d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool);
274d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*);
275014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool );
276f746d493SDmitry Karpeev 
2773d957683SBarry Smith /*E
2783d957683SBarry Smith     PCASMType - Type of additive Schwarz method to use
2793d957683SBarry Smith 
280feb221f8SDmitry Karpeev $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
281feb221f8SDmitry Karpeev $                        and computed values in ghost regions are added together.
282feb221f8SDmitry Karpeev $                        Classical standard additive Schwarz.
283feb221f8SDmitry Karpeev $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
284feb221f8SDmitry Karpeev $                        region are discarded.
285feb221f8SDmitry Karpeev $                        Default.
286feb221f8SDmitry Karpeev $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
287feb221f8SDmitry Karpeev $                        region are added back in.
288feb221f8SDmitry Karpeev $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
289feb221f8SDmitry Karpeev $                        discarded.
290feb221f8SDmitry Karpeev $                        Not very good.
2913d957683SBarry Smith 
2923d957683SBarry Smith    Level: beginner
2933d957683SBarry Smith 
2943d957683SBarry Smith .seealso: PCASMSetType()
2953d957683SBarry Smith E*/
296d252947aSBarry Smith typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
2976a6fc655SJed Brown PETSC_EXTERN const char *const PCASMTypes[];
2989dcbbd2bSBarry Smith 
299014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType);
300014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]);
301014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]);
302014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
303014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]);
304014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]);
305981c4779SBarry Smith 
3063d957683SBarry Smith /*E
3076a4f0f73SDmitry Karpeev     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
308f746d493SDmitry Karpeev 
3096a4f0f73SDmitry Karpeev    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
3106a4f0f73SDmitry Karpeev    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
3116a4f0f73SDmitry Karpeev    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
3126a4f0f73SDmitry Karpeev    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
3136a4f0f73SDmitry Karpeev 
3146a4f0f73SDmitry Karpeev $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
3156a4f0f73SDmitry Karpeev $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
3166a4f0f73SDmitry Karpeev $                        from neighboring subdomains.
317feb221f8SDmitry Karpeev $                        Classical standard additive Schwarz.
3186a4f0f73SDmitry Karpeev $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
3196a4f0f73SDmitry Karpeev $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
3206a4f0f73SDmitry Karpeev $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
3216a4f0f73SDmitry Karpeev $                        assumption).
322feb221f8SDmitry Karpeev $                        Default.
3236a4f0f73SDmitry Karpeev $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
3246a4f0f73SDmitry Karpeev $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
3256a4f0f73SDmitry Karpeev $                        from neighboring subdomains.
3266a4f0f73SDmitry Karpeev $
3276a4f0f73SDmitry Karpeev $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
328feb221f8SDmitry Karpeev $                        Not very good.
329f746d493SDmitry Karpeev 
330f746d493SDmitry Karpeev    Level: beginner
331f746d493SDmitry Karpeev 
332f746d493SDmitry Karpeev .seealso: PCGASMSetType()
333f746d493SDmitry Karpeev E*/
334f746d493SDmitry Karpeev typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
3356a6fc655SJed Brown PETSC_EXTERN const char *const PCGASMTypes[];
336f746d493SDmitry Karpeev 
337014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]);
338014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt,PetscBool);
339014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt);
340d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCGASMSetDMSubdomains(PC,PetscBool);
341d709ea83SDmitry Karpeev PETSC_EXTERN PetscErrorCode PCGASMGetDMSubdomains(PC,PetscBool*);
342014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool );
343f746d493SDmitry Karpeev 
344014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType);
345014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]);
346014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]);
347014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**);
348014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]);
349014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]);
350f746d493SDmitry Karpeev 
351f746d493SDmitry Karpeev /*E
3523d957683SBarry Smith     PCCompositeType - Determines how two or more preconditioner are composed
3533d957683SBarry Smith 
3543d957683SBarry Smith $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
3553d957683SBarry Smith $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
3563d957683SBarry Smith $                                computed after the previous preconditioner application
357ccb205f8SBarry Smith $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
358ccb205f8SBarry Smith $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
3593d957683SBarry Smith $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
3603d957683SBarry Smith $                         where first preconditioner is built from alpha I + S and second from
3613d957683SBarry Smith $                         alpha I + R
3623d957683SBarry Smith 
3633d957683SBarry Smith    Level: beginner
3643d957683SBarry Smith 
3653d957683SBarry Smith .seealso: PCCompositeSetType()
3663d957683SBarry Smith E*/
3673b224e63SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
3686a6fc655SJed Brown PETSC_EXTERN const char *const PCCompositeTypes[];
3699dcbbd2bSBarry Smith 
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType);
37119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar);
374981c4779SBarry Smith 
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt);
376014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter);
377014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*);
378da3a660dSBarry Smith 
379014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double);
380014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt);
381014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt);
382014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt);
383014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt);
384014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt);
385014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt);
386014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt);
3873304466cSBarry Smith 
388014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]);
389014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]);
390014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]);
391014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]);
3923304466cSBarry Smith 
393014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*);
394014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*);
395014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType);
396014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt);
397014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS);
398014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*);
3994ab8060aSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool);
4004ab8060aSDmitry Karpeev PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*);
401084e4875SJed Brown 
402084e4875SJed Brown /*E
403084e4875SJed Brown     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
404084e4875SJed Brown 
405084e4875SJed Brown     Level: intermediate
406084e4875SJed Brown 
407084e4875SJed Brown .seealso: PCFieldSplitSchurPrecondition()
408084e4875SJed Brown E*/
409e87fab1bSBarry Smith typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_A11,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType;
410014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
411084e4875SJed Brown 
412ab1df9f5SJed Brown /*E
413c9c6ffaaSJed Brown     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
414ab1df9f5SJed Brown 
415ab1df9f5SJed Brown     Level: intermediate
416ab1df9f5SJed Brown 
417c9c6ffaaSJed Brown .seealso: PCFieldSplitSetSchurFactType()
418ab1df9f5SJed Brown E*/
419ab1df9f5SJed Brown typedef enum {
420c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
421c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
422c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
423c9c6ffaaSJed Brown   PC_FIELDSPLIT_SCHUR_FACT_FULL
424c9c6ffaaSJed Brown } PCFieldSplitSchurFactType;
425014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
426ab1df9f5SJed Brown 
427014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat);
428014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType);
429014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*);
4303d30b1ffSBarry Smith 
431014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat);
432014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat);
4332a6744ebSBarry Smith 
434014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*);
435014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *);
4362102ba4dSSatish Balay 
437014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]);
43867fac13cSBarry Smith 
439014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM);
440014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*);
4416c699258SBarry Smith 
442014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*);
443014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*);
4441b2093e4SBarry Smith 
445014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal);
446014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt);
447014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool);
4486c699258SBarry Smith 
449014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal);
450014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool);
451014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt);
452014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt);
45337f80224SJose E. Roman /*E
45437f80224SJose E. Roman     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
45537f80224SJose E. Roman 
45637f80224SJose E. Roman     Level: intermediate
45737f80224SJose E. Roman 
45837f80224SJose E. Roman .seealso: PCPARMSSetGlobal()
45937f80224SJose E. Roman E*/
46037f80224SJose E. Roman typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
4616a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
46237f80224SJose E. Roman /*E
46337f80224SJose E. Roman     PCPARMSLocalType - Determines the local preconditioner method in PARMS
46437f80224SJose E. Roman 
46537f80224SJose E. Roman     Level: intermediate
46637f80224SJose E. Roman 
46737f80224SJose E. Roman .seealso: PCPARMSSetLocal()
46837f80224SJose E. Roman E*/
46937f80224SJose E. Roman typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
4706a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSLocalTypes[];
47137f80224SJose E. Roman 
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type);
473014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type);
474014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits);
475014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart);
476014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym);
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2);
47837f80224SJose E. Roman 
479014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt);
480014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool);
481014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool);
482014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt);
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal);
484014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt);
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt);
48619fd82e9SBarry Smith typedef const char* PCGAMGType;
48719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType );
488014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n);
489014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n);
490014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool);
491fc897844SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool);
49256de70efSSatish Balay 
4930c7d97c5SJed Brown #if defined(PETSC_HAVE_PCBDDC)
49453cdbc3dSStefano Zampini /* Enum defining how to treat the coarse problem */
4950c7d97c5SJed Brown typedef enum {SEQUENTIAL_BDDC,REPLICATED_BDDC,PARALLEL_BDDC,MULTILEVEL_BDDC} CoarseProblemType;
4964fad6a16SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt);
4974fad6a16SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetMaxLevels(PC,PetscInt);
4980bdf917eSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace);
499014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS);
500014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*);
501014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS);
502014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*);
503014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseProblemType(PC,CoarseProblemType);
504014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]);
5051a83f524SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode);
5063425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*);
5073425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec);
5083425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec);
5090c7d97c5SJed Brown #endif
5100c7d97c5SJed Brown 
511b965d45eSStefano Zampini PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool);
512014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar);
51346caae47SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec);
514d3b1e0e7SMatthew Knepley 
515b4876bcbSBarry Smith /*E
516b4876bcbSBarry Smith     PCMGType - Determines the type of multigrid method that is run.
517b4876bcbSBarry Smith 
518b4876bcbSBarry Smith    Level: beginner
519b4876bcbSBarry Smith 
520b4876bcbSBarry Smith    Values:
521b4876bcbSBarry Smith +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles()
522b4876bcbSBarry Smith .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
523b4876bcbSBarry Smith                 smoothed before updating the residual. This only uses the
524b4876bcbSBarry Smith                 down smoother, in the preconditioner the upper smoother is ignored
525b4876bcbSBarry Smith .  PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
526b4876bcbSBarry Smith             that is starts on the coarsest grid, performs a cycle, interpolates
527b4876bcbSBarry Smith             to the next, performs a cycle etc. This is much like the F-cycle presented in "Multigrid" by Trottenberg, Oosterlee, Schuller page 49, but that
528b4876bcbSBarry Smith             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
529b4876bcbSBarry Smith             calls the V-cycle only on the coarser level and has a post-smoother instead.
530b4876bcbSBarry Smith -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
531b4876bcbSBarry Smith                from a finer
532b4876bcbSBarry Smith 
533b4876bcbSBarry Smith .seealso: PCMGSetType()
534b4876bcbSBarry Smith 
535b4876bcbSBarry Smith E*/
536b4876bcbSBarry Smith typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
537b4876bcbSBarry Smith PETSC_EXTERN const char *const PCMGTypes[];
538b4876bcbSBarry Smith #define PC_MG_CASCADE PC_MG_KASKADE;
539b4876bcbSBarry Smith 
540b4876bcbSBarry Smith /*E
541b4876bcbSBarry Smith     PCMGCycleType - Use V-cycle or W-cycle
542b4876bcbSBarry Smith 
543b4876bcbSBarry Smith    Level: beginner
544b4876bcbSBarry Smith 
545b4876bcbSBarry Smith    Values:
546b4876bcbSBarry Smith +  PC_MG_V_CYCLE
547b4876bcbSBarry Smith -  PC_MG_W_CYCLE
548b4876bcbSBarry Smith 
549b4876bcbSBarry Smith .seealso: PCMGSetCycleType()
550b4876bcbSBarry Smith 
551b4876bcbSBarry Smith E*/
552b4876bcbSBarry Smith typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
553b4876bcbSBarry Smith PETSC_EXTERN const char *const PCMGCycleTypes[];
554b4876bcbSBarry Smith 
555b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType);
556b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*);
557b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*);
558b4876bcbSBarry Smith 
559b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt);
560b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt);
561b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType);
562b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType);
563b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt);
564b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt);
565b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool);
566b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*);
567b4876bcbSBarry Smith 
568b4876bcbSBarry Smith 
569b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec);
570b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec);
571b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec);
572b4876bcbSBarry Smith 
573b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat);
574b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*);
575b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat);
576b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*);
577b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec);
578b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*);
579b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat);
580b4876bcbSBarry Smith 
581b4876bcbSBarry Smith /*E
582b4876bcbSBarry Smith     PCExoticType - Face based or wirebasket based coarse grid space
583b4876bcbSBarry Smith 
584b4876bcbSBarry Smith    Level: beginner
585b4876bcbSBarry Smith 
586b4876bcbSBarry Smith .seealso: PCExoticSetType(), PCEXOTIC
587b4876bcbSBarry Smith E*/
588b4876bcbSBarry Smith typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
589b4876bcbSBarry Smith PETSC_EXTERN const char *const PCExoticTypes[];
590b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);
591b4876bcbSBarry Smith 
592123ea438SMatthew Knepley #endif /* __PETSCPC_H */
593