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 */ 14140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList 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 11155849f57SBarry Smith #define PC_FILE_CLASSID 1211222 11255849f57SBarry Smith 1134d0a8057SBarry Smith /*E 1144d0a8057SBarry Smith PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates 1154d0a8057SBarry Smith 1164d0a8057SBarry Smith Level: advanced 1174d0a8057SBarry Smith 1184d0a8057SBarry Smith Notes: this must match finclude/petscpc.h and the KSPConvergedReason values in petscksp.h 1194d0a8057SBarry Smith 1204d0a8057SBarry Smith .seealso: PCApplyRichardson() 1214d0a8057SBarry Smith E*/ 1224d0a8057SBarry Smith typedef enum { 1234d0a8057SBarry Smith PCRICHARDSON_CONVERGED_RTOL = 2, 1244d0a8057SBarry Smith PCRICHARDSON_CONVERGED_ATOL = 3, 1254d0a8057SBarry Smith PCRICHARDSON_CONVERGED_ITS = 4, 1264d0a8057SBarry Smith PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason; 1274d0a8057SBarry Smith 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*); 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *); 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetInitialGuessNonzero(PC,PetscBool ); 13184cb2905SBarry Smith 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRegisterDestroy(void); 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRegisterAll(const char[]); 134014dd563SJed Brown PETSC_EXTERN PetscBool PCRegisterAllCalled; 13584cb2905SBarry Smith 136014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRegister(const char[],const char[],const char[],PetscErrorCode(*)(PC)); 13730de9b25SBarry Smith 13830de9b25SBarry Smith /*MC 13930de9b25SBarry Smith PCRegisterDynamic - Adds a method to the preconditioner package. 14030de9b25SBarry Smith 14130de9b25SBarry Smith Synopsis: 142f2ba6396SBarry Smith #include "petscpc.h" 1431890ba74SBarry Smith PetscErrorCode PCRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(PC)) 14430de9b25SBarry Smith 14530de9b25SBarry Smith Not collective 14630de9b25SBarry Smith 14730de9b25SBarry Smith Input Parameters: 14830de9b25SBarry Smith + name_solver - name of a new user-defined solver 14930de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 15030de9b25SBarry Smith . name_create - name of routine to create method context 15130de9b25SBarry Smith - routine_create - routine to create method context 15230de9b25SBarry Smith 15330de9b25SBarry Smith Notes: 15430de9b25SBarry Smith PCRegisterDynamic() may be called multiple times to add several user-defined preconditioners. 15530de9b25SBarry Smith 15630de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 15730de9b25SBarry Smith is ignored. 15830de9b25SBarry Smith 15930de9b25SBarry Smith Sample usage: 16030de9b25SBarry Smith .vb 16130de9b25SBarry Smith PCRegisterDynamic("my_solver","/home/username/my_lib/lib/libO/solaris/mylib", 16230de9b25SBarry Smith "MySolverCreate",MySolverCreate); 16330de9b25SBarry Smith .ve 16430de9b25SBarry Smith 16530de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 16630de9b25SBarry Smith $ PCSetType(pc,"my_solver") 16730de9b25SBarry Smith or at runtime via the option 16830de9b25SBarry Smith $ -pc_type my_solver 16930de9b25SBarry Smith 17030de9b25SBarry Smith Level: advanced 17130de9b25SBarry Smith 172ab901514SBarry Smith Notes: ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, or ${any environmental variable} 17330de9b25SBarry Smith occuring in pathname will be replaced with appropriate values. 17430de9b25SBarry Smith If your function is not being put into a shared library then use PCRegister() instead 17530de9b25SBarry Smith 17630de9b25SBarry Smith .keywords: PC, register 17730de9b25SBarry Smith 17830de9b25SBarry Smith .seealso: PCRegisterAll(), PCRegisterDestroy() 17930de9b25SBarry Smith M*/ 180aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 181f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,0) 1826df38c32SLois Curfman McInnes #else 183f1af5d2fSBarry Smith #define PCRegisterDynamic(a,b,c,d) PCRegister(a,b,c,d) 1846df38c32SLois Curfman McInnes #endif 1856df38c32SLois Curfman McInnes 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCReset(PC); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDestroy(PC*); 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC); 18919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*); 19014c91fddSBarry Smith 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*); 192014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*); 193014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*); 1945b116368SBarry Smith 195014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat,MatStructure); 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*,MatStructure*); 197014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *); 1984b0e389bSBarry Smith 199014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer); 20055849f57SBarry Smith PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer); 2017bc3d0afSSatish Balay 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]); 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]); 2058ed539a5SBarry Smith 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*); 20771601f6fSBarry Smith 208d6913704SBarry Smith /* 209d6913704SBarry Smith These are used to provide extra scaling of preconditioned 2100f3b3ca1SHong Zhang operator for time-stepping schemes like in SUNDIALS 211d6913704SBarry Smith */ 212014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *); 213014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec); 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec); 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec); 216d6913704SBarry Smith 21784cb2905SBarry Smith /* ------------- options specific to particular preconditioners --------- */ 218329f5518SBarry Smith 219014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowMax(PC); 220014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseRowSum(PC); 221014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC); 222014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType); 223014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal); 224014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt); 225d03aef70SBarry Smith 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal); 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCEisenstatNoDiagonalScaling(PC); 228421c37bdSBarry Smith 22915aa81f8SBarry Smith #define USE_PRECONDITIONER_MATRIX 0 23015aa81f8SBarry Smith #define USE_TRUE_MATRIX 1 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiSetUseTrueLocal(PC); 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 PCKSPSetUseTrue(PC); 236981c4779SBarry Smith 237014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec)); 238014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)); 239014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec)); 240014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC)); 241014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*)); 242014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer)); 243014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC)); 244014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**); 245014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*); 246014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]); 247014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]); 248aabeff55SBarry Smith 249014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal); 250d90ac83dSHong Zhang 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType); 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal); 253d90ac83dSHong Zhang 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverPackage(PC,const MatSolverPackage); 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverPackage(PC,const MatSolverPackage*); 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverPackage(PC); 2572401956bSBarry Smith 258014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal); 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal); 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal); 261421c37bdSBarry Smith 26219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool ); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool ); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC); 266014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC); 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool ); 268f5a88c2aSLois Curfman McInnes 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt); 271b35a507dSBarry Smith 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]); 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]); 274014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt); 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); 340014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool ); 341f746d493SDmitry Karpeev 342014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType); 343014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat,PetscInt,PetscInt,IS*[],IS*[]); 344014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS[],IS[]); 345014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 346014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]); 347014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]); 348f746d493SDmitry Karpeev 349f746d493SDmitry Karpeev /*E 3503d957683SBarry Smith PCCompositeType - Determines how two or more preconditioner are composed 3513d957683SBarry Smith 3523d957683SBarry Smith $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 3533d957683SBarry Smith $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 3543d957683SBarry Smith $ computed after the previous preconditioner application 355ccb205f8SBarry Smith $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 356ccb205f8SBarry Smith $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions) 3573d957683SBarry Smith $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 3583d957683SBarry Smith $ where first preconditioner is built from alpha I + S and second from 3593d957683SBarry Smith $ alpha I + R 3603d957683SBarry Smith 3613d957683SBarry Smith Level: beginner 3623d957683SBarry Smith 3633d957683SBarry Smith .seealso: PCCompositeSetType() 3643d957683SBarry Smith E*/ 3653b224e63SBarry Smith typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType; 3666a6fc655SJed Brown PETSC_EXTERN const char *const PCCompositeTypes[]; 3679dcbbd2bSBarry Smith 368014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSetUseTrue(PC); 369014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType); 37019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType); 371014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *); 372014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar); 373981c4779SBarry Smith 374014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt); 375014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter); 376014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*); 377da3a660dSBarry Smith 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double); 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt); 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt); 3863304466cSBarry Smith 387014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]); 388014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]); 389014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]); 390014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]); 3913304466cSBarry Smith 392014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*); 393014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*); 394014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType); 395014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt); 396014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS); 397014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*); 398084e4875SJed Brown 399084e4875SJed Brown /*E 400084e4875SJed Brown PCFieldSplitSchurPreType - Determines how to precondition Schur complement 401084e4875SJed Brown 402084e4875SJed Brown Level: intermediate 403084e4875SJed Brown 404084e4875SJed Brown .seealso: PCFieldSplitSchurPrecondition() 405084e4875SJed Brown E*/ 406e87fab1bSBarry Smith typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_A11,PC_FIELDSPLIT_SCHUR_PRE_USER} PCFieldSplitSchurPreType; 407014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[]; 408084e4875SJed Brown 409ab1df9f5SJed Brown /*E 410c9c6ffaaSJed Brown PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 411ab1df9f5SJed Brown 412ab1df9f5SJed Brown Level: intermediate 413ab1df9f5SJed Brown 414c9c6ffaaSJed Brown .seealso: PCFieldSplitSetSchurFactType() 415ab1df9f5SJed Brown E*/ 416ab1df9f5SJed Brown typedef enum { 417c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_DIAG, 418c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_LOWER, 419c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_UPPER, 420c9c6ffaaSJed Brown PC_FIELDSPLIT_SCHUR_FACT_FULL 421c9c6ffaaSJed Brown } PCFieldSplitSchurFactType; 422014dd563SJed Brown PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[]; 423ab1df9f5SJed Brown 424014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat); 425014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType); 426014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*); 4273d30b1ffSBarry Smith 428014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat); 429014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat); 4302a6744ebSBarry Smith 431014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*); 432014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSASetVectors(PC,PetscInt,PetscReal *); 4332102ba4dSSatish Balay 434014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]); 43567fac13cSBarry Smith 436014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM); 437014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*); 4386c699258SBarry Smith 439014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*); 440014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*); 4411b2093e4SBarry Smith 442014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetTolerance(PC,PetscReal); 443014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetIterations(PC,PetscInt); 444014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBiCGStabCUSPSetUseVerboseMonitor(PC,PetscBool); 4456c699258SBarry Smith 446014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetDropTolerance(PC,PetscReal); 447014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPUseScaling(PC,PetscBool); 448014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetNonzeros(PC,PetscInt); 449014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCAINVCUSPSetLinParameter(PC,PetscInt); 45037f80224SJose E. Roman /*E 45137f80224SJose E. Roman PCPARMSGlobalType - Determines the global preconditioner method in PARMS 45237f80224SJose E. Roman 45337f80224SJose E. Roman Level: intermediate 45437f80224SJose E. Roman 45537f80224SJose E. Roman .seealso: PCPARMSSetGlobal() 45637f80224SJose E. Roman E*/ 45737f80224SJose E. Roman typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType; 4586a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSGlobalTypes[]; 45937f80224SJose E. Roman /*E 46037f80224SJose E. Roman PCPARMSLocalType - Determines the local preconditioner method in PARMS 46137f80224SJose E. Roman 46237f80224SJose E. Roman Level: intermediate 46337f80224SJose E. Roman 46437f80224SJose E. Roman .seealso: PCPARMSSetLocal() 46537f80224SJose E. Roman E*/ 46637f80224SJose E. Roman typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType; 4676a6fc655SJed Brown PETSC_EXTERN const char *const PCPARMSLocalTypes[]; 46837f80224SJose E. Roman 469014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type); 470014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type); 471014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits); 472014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart); 473014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym); 474014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2); 47537f80224SJose E. Roman 476014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt); 477014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetRepartitioning(PC,PetscBool); 478014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetUseASMAggs(PC,PetscBool); 479014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt); 480014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal); 481014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt); 482014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt); 48319fd82e9SBarry Smith typedef const char* PCGAMGType; 48419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType ); 485014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC pc, PetscInt n); 486014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC pc, PetscBool n); 487014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscBool); 488fc897844SJed Brown PETSC_EXTERN PetscErrorCode PCGAMGSetReuseProl(PC,PetscBool); 48956de70efSSatish Balay 4900c7d97c5SJed Brown #if defined(PETSC_HAVE_PCBDDC) 49153cdbc3dSStefano Zampini /* Enum defining how to treat the coarse problem */ 4920c7d97c5SJed Brown typedef enum {SEQUENTIAL_BDDC,REPLICATED_BDDC,PARALLEL_BDDC,MULTILEVEL_BDDC} CoarseProblemType; 4934fad6a16SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt); 4944fad6a16SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetMaxLevels(PC,PetscInt); 4950bdf917eSStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCSetNullSpace(PC,MatNullSpace); 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS); 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*); 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS); 499014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*); 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseProblemType(PC,CoarseProblemType); 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]); 5021a83f524SJed Brown PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode); 5033425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,Mat*,PC*); 5043425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec); 5053425bc38SStefano Zampini PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec); 5060c7d97c5SJed Brown #endif 5070c7d97c5SJed Brown 508b965d45eSStefano Zampini PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool); 509014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar); 51046caae47SJed Brown PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec); 511d3b1e0e7SMatthew Knepley 512*b4876bcbSBarry Smith /*E 513*b4876bcbSBarry Smith PCMGType - Determines the type of multigrid method that is run. 514*b4876bcbSBarry Smith 515*b4876bcbSBarry Smith Level: beginner 516*b4876bcbSBarry Smith 517*b4876bcbSBarry Smith Values: 518*b4876bcbSBarry Smith + PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles() 519*b4876bcbSBarry Smith . PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are 520*b4876bcbSBarry Smith smoothed before updating the residual. This only uses the 521*b4876bcbSBarry Smith down smoother, in the preconditioner the upper smoother is ignored 522*b4876bcbSBarry Smith . PC_MG_FULL - same as multiplicative except one also performs grid sequencing, 523*b4876bcbSBarry Smith that is starts on the coarsest grid, performs a cycle, interpolates 524*b4876bcbSBarry 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 525*b4876bcbSBarry Smith algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 526*b4876bcbSBarry Smith calls the V-cycle only on the coarser level and has a post-smoother instead. 527*b4876bcbSBarry Smith - PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level 528*b4876bcbSBarry Smith from a finer 529*b4876bcbSBarry Smith 530*b4876bcbSBarry Smith .seealso: PCMGSetType() 531*b4876bcbSBarry Smith 532*b4876bcbSBarry Smith E*/ 533*b4876bcbSBarry Smith typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType; 534*b4876bcbSBarry Smith PETSC_EXTERN const char *const PCMGTypes[]; 535*b4876bcbSBarry Smith #define PC_MG_CASCADE PC_MG_KASKADE; 536*b4876bcbSBarry Smith 537*b4876bcbSBarry Smith /*E 538*b4876bcbSBarry Smith PCMGCycleType - Use V-cycle or W-cycle 539*b4876bcbSBarry Smith 540*b4876bcbSBarry Smith Level: beginner 541*b4876bcbSBarry Smith 542*b4876bcbSBarry Smith Values: 543*b4876bcbSBarry Smith + PC_MG_V_CYCLE 544*b4876bcbSBarry Smith - PC_MG_W_CYCLE 545*b4876bcbSBarry Smith 546*b4876bcbSBarry Smith .seealso: PCMGSetCycleType() 547*b4876bcbSBarry Smith 548*b4876bcbSBarry Smith E*/ 549*b4876bcbSBarry Smith typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType; 550*b4876bcbSBarry Smith PETSC_EXTERN const char *const PCMGCycleTypes[]; 551*b4876bcbSBarry Smith 552*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType); 553*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*); 554*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*); 555*b4876bcbSBarry Smith 556*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothUp(PC,PetscInt); 557*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetNumberSmoothDown(PC,PetscInt); 558*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType); 559*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType); 560*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetCyclesOnLevel(PC,PetscInt,PetscInt); 561*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt); 562*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PetscBool); 563*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PetscBool*); 564*b4876bcbSBarry Smith 565*b4876bcbSBarry Smith 566*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec); 567*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec); 568*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec); 569*b4876bcbSBarry Smith 570*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat); 571*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*); 572*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat); 573*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*); 574*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec); 575*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*); 576*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat); 577*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGDefaultResidual(Mat,Vec,Vec,Vec); 578*b4876bcbSBarry Smith 579*b4876bcbSBarry Smith /*E 580*b4876bcbSBarry Smith PCExoticType - Face based or wirebasket based coarse grid space 581*b4876bcbSBarry Smith 582*b4876bcbSBarry Smith Level: beginner 583*b4876bcbSBarry Smith 584*b4876bcbSBarry Smith .seealso: PCExoticSetType(), PCEXOTIC 585*b4876bcbSBarry Smith E*/ 586*b4876bcbSBarry Smith typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType; 587*b4876bcbSBarry Smith PETSC_EXTERN const char *const PCExoticTypes[]; 588*b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType); 589*b4876bcbSBarry Smith 590123ea438SMatthew Knepley #endif /* __PETSCPC_H */ 591