1 /* 2 Preconditioner module. 3 */ 4 #if !defined(__PETSCPC_H) 5 #define __PETSCPC_H 6 #include <petscmat.h> 7 #include <petscpctypes.h> 8 9 PETSC_EXTERN PetscErrorCode PCInitializePackage(void); 10 11 /* 12 PCList contains the list of preconditioners currently registered 13 These are added with PCRegister() 14 */ 15 PETSC_EXTERN PetscFunctionList PCList; 16 17 /* Logging support */ 18 PETSC_EXTERN PetscClassId PC_CLASSID; 19 20 PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm,PC*); 21 PETSC_EXTERN PetscErrorCode PCSetType(PC,PCType); 22 PETSC_EXTERN PetscErrorCode PCGetType(PC,PCType*); 23 PETSC_EXTERN PetscErrorCode PCSetUp(PC); 24 PETSC_EXTERN PetscErrorCode PCGetSetUpFailedReason(PC,PCFailedReason*); 25 PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC); 26 PETSC_EXTERN PetscErrorCode PCApply(PC,Vec,Vec); 27 PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC,Vec,Vec); 28 PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC,Vec,Vec); 29 PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC,PCSide,Vec,Vec,Vec); 30 PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC,Vec,Vec); 31 PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC,PetscBool *); 32 PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC,PCSide,Vec,Vec,Vec); 33 PETSC_EXTERN PetscErrorCode PCSetReusePreconditioner(PC,PetscBool); 34 PETSC_EXTERN PetscErrorCode PCGetReusePreconditioner(PC,PetscBool*); 35 PETSC_EXTERN PetscErrorCode PCSetErrorIfFailure(PC,PetscBool); 36 37 #define PC_FILE_CLASSID 1211222 38 39 PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*); 40 PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC,PetscBool *); 41 PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC,PetscBool); 42 PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC,PetscBool*); 43 44 45 PETSC_EXTERN PetscErrorCode PCRegister(const char[],PetscErrorCode(*)(PC)); 46 47 PETSC_EXTERN PetscErrorCode PCReset(PC); 48 PETSC_EXTERN PetscErrorCode PCDestroy(PC*); 49 PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC); 50 51 PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC,Mat*); 52 PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC,PetscErrorCode(*)(PC,PetscInt,const IS[],const IS[],Mat[],void*),void*); 53 PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC,PetscInt,const IS[],const IS[],Mat[],void*); 54 55 PETSC_EXTERN PetscErrorCode PCSetOperators(PC,Mat,Mat); 56 PETSC_EXTERN PetscErrorCode PCGetOperators(PC,Mat*,Mat*); 57 PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC,PetscBool *,PetscBool *); 58 59 PETSC_EXTERN PetscErrorCode PCView(PC,PetscViewer); 60 PETSC_EXTERN PetscErrorCode PCLoad(PC,PetscViewer); 61 PETSC_STATIC_INLINE PetscErrorCode PCViewFromOptions(PC A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 62 63 PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC,const char[]); 64 PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC,const char[]); 65 PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC,const char*[]); 66 67 PETSC_EXTERN PetscErrorCode PCComputeExplicitOperator(PC,Mat*); 68 69 /* 70 These are used to provide extra scaling of preconditioned 71 operator for time-stepping schemes like in SUNDIALS 72 */ 73 PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC,PetscBool *); 74 PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC,Vec,Vec); 75 PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC,Vec,Vec); 76 PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC,Vec); 77 78 /* ------------- options specific to particular preconditioners --------- */ 79 80 PETSC_EXTERN PetscErrorCode PCJacobiSetType(PC,PCJacobiType); 81 PETSC_EXTERN PetscErrorCode PCJacobiGetType(PC,PCJacobiType*); 82 PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC,PetscBool); 83 PETSC_EXTERN PetscErrorCode PCJacobiGetUseAbs(PC,PetscBool*); 84 PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC,MatSORType); 85 PETSC_EXTERN PetscErrorCode PCSORGetSymmetric(PC,MatSORType*); 86 PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC,PetscReal); 87 PETSC_EXTERN PetscErrorCode PCSORGetOmega(PC,PetscReal*); 88 PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC,PetscInt,PetscInt); 89 PETSC_EXTERN PetscErrorCode PCSORGetIterations(PC,PetscInt*,PetscInt*); 90 91 PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC,PetscReal); 92 PETSC_EXTERN PetscErrorCode PCEisenstatGetOmega(PC,PetscReal*); 93 PETSC_EXTERN PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC,PetscBool); 94 PETSC_EXTERN PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC,PetscBool*); 95 96 PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC,PetscInt,const PetscInt[]); 97 PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC,PetscInt*,const PetscInt*[]); 98 PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC,PetscInt,const PetscInt[]); 99 PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC,PetscInt*,const PetscInt*[]); 100 101 PETSC_EXTERN PetscErrorCode PCShellSetApply(PC,PetscErrorCode (*)(PC,Vec,Vec)); 102 PETSC_EXTERN PetscErrorCode PCShellSetApplySymmetricLeft(PC,PetscErrorCode (*)(PC,Vec,Vec)); 103 PETSC_EXTERN PetscErrorCode PCShellSetApplySymmetricRight(PC,PetscErrorCode (*)(PC,Vec,Vec)); 104 PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC,PetscErrorCode (*)(PC,PCSide,Vec,Vec,Vec)); 105 PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC,PetscErrorCode (*)(PC,Vec,Vec)); 106 PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC,PetscErrorCode (*)(PC)); 107 PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC,PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,PetscInt,PetscBool ,PetscInt*,PCRichardsonConvergedReason*)); 108 PETSC_EXTERN PetscErrorCode PCShellSetView(PC,PetscErrorCode (*)(PC,PetscViewer)); 109 PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC,PetscErrorCode (*)(PC)); 110 PETSC_EXTERN PetscErrorCode PCShellSetContext(PC,void*); 111 PETSC_EXTERN PetscErrorCode PCShellGetContext(PC,void**); 112 PETSC_EXTERN PetscErrorCode PCShellSetName(PC,const char[]); 113 PETSC_EXTERN PetscErrorCode PCShellGetName(PC,const char*[]); 114 115 PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC,PetscReal); 116 117 PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC,MatFactorShiftType); 118 PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC,PetscReal); 119 120 PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverType(PC,MatSolverType); 121 PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverType(PC,MatSolverType*); 122 PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverType(PC); 123 124 PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC,PetscReal); 125 PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC,PetscReal); 126 PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC,PetscReal); 127 128 PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC,MatOrderingType); 129 PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC,PetscBool ); 130 PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC,PetscBool ); 131 PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC,PetscBool); 132 PETSC_EXTERN PetscErrorCode PCFactorGetUseInPlace(PC,PetscBool*); 133 PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC,PetscBool); 134 PETSC_EXTERN PetscErrorCode PCFactorGetAllowDiagonalFill(PC,PetscBool*); 135 PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC,PetscBool); 136 137 PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC,PetscInt); 138 PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC,PetscInt*); 139 PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC,PetscReal,PetscReal,PetscInt); 140 PETSC_EXTERN PetscErrorCode PCFactorGetZeroPivot(PC,PetscReal*); 141 PETSC_EXTERN PetscErrorCode PCFactorGetShiftAmount(PC,PetscReal*); 142 PETSC_EXTERN PetscErrorCode PCFactorGetShiftType(PC,MatFactorShiftType*); 143 144 PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC,PetscInt,IS[],IS[]); 145 PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC,PetscInt,IS[],IS[]); 146 PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC,PetscInt); 147 PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC,PetscBool); 148 PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC,PetscBool*); 149 PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC,PetscBool); 150 151 PETSC_EXTERN PetscErrorCode PCASMSetType(PC,PCASMType); 152 PETSC_EXTERN PetscErrorCode PCASMGetType(PC,PCASMType*); 153 PETSC_EXTERN PetscErrorCode PCASMSetLocalType(PC,PCCompositeType); 154 PETSC_EXTERN PetscErrorCode PCASMGetLocalType(PC,PCCompositeType*); 155 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat,PetscInt,IS*[]); 156 PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt,IS[],IS[]); 157 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 158 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC,PetscInt*,IS*[],IS*[]); 159 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC,PetscInt*,Mat*[]); 160 PETSC_EXTERN PetscErrorCode PCASMGetSubMatType(PC,MatType*); 161 PETSC_EXTERN PetscErrorCode PCASMSetSubMatType(PC,MatType); 162 163 PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC,PetscInt); 164 PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC,PetscInt,IS[],IS[]); 165 PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC,PetscInt); 166 PETSC_EXTERN PetscErrorCode PCGASMSetUseDMSubdomains(PC,PetscBool); 167 PETSC_EXTERN PetscErrorCode PCGASMGetUseDMSubdomains(PC,PetscBool*); 168 PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC,PetscBool ); 169 170 PETSC_EXTERN PetscErrorCode PCGASMSetType(PC,PCGASMType); 171 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains(Mat,PetscInt,PetscInt*,IS*[]); 172 PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt,IS*[],IS*[]); 173 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,IS**,IS**); 174 PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC,PetscInt*,IS*[],IS*[]); 175 PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC,PetscInt*,Mat*[]); 176 177 PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC,PCCompositeType); 178 PETSC_EXTERN PetscErrorCode PCCompositeGetType(PC,PCCompositeType*); 179 PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC,PCType); 180 PETSC_EXTERN PetscErrorCode PCCompositeGetNumberPC(PC,PetscInt *); 181 PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC,PetscInt,PC *); 182 PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC,PetscScalar); 183 184 PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC,PetscInt); 185 PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC,VecScatter,VecScatter); 186 PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC,Mat*,Mat*); 187 188 PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC,double); 189 PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC,PetscInt); 190 PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC,PetscInt); 191 PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC,PetscInt); 192 PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC,PetscInt); 193 PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC,PetscInt); 194 PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC,PetscInt); 195 PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC,PetscInt); 196 197 PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC,const char[]); 198 PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC,const char*[]); 199 PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteGradient(PC,Mat); 200 PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteCurl(PC,Mat); 201 PETSC_EXTERN PetscErrorCode PCHYPRESetInterpolations(PC,PetscInt,Mat,Mat[],Mat,Mat[]); 202 PETSC_EXTERN PetscErrorCode PCHYPRESetEdgeConstantVectors(PC,Vec,Vec,Vec); 203 PETSC_EXTERN PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC,Mat); 204 PETSC_EXTERN PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC,Mat); 205 206 PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC,const char[],PetscInt,const PetscInt*,const PetscInt*); 207 PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC,PCCompositeType); 208 PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC,PCCompositeType*); 209 PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC,PetscInt); 210 PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC,const char[],IS); 211 PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC,const char[],IS*); 212 PETSC_EXTERN PetscErrorCode PCFieldSplitRestrictIS(PC,IS); 213 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC,PetscBool); 214 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC,PetscBool*); 215 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC,PetscBool); 216 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC,PetscBool*); 217 PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC,PetscBool); 218 PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC,PetscBool*); 219 220 PETSC_EXTERN PETSC_DEPRECATED("Use PCFieldSplitSetSchurPre") PetscErrorCode PCFieldSplitSchurPrecondition(PC,PCFieldSplitSchurPreType,Mat); 221 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC,PCFieldSplitSchurPreType,Mat); 222 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC,PCFieldSplitSchurPreType*,Mat*); 223 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC,PCFieldSplitSchurFactType); 224 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurScale(PC,PetscScalar); 225 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC,Mat*,Mat*,Mat*,Mat*); 226 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC,Mat *S); 227 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC,Mat *S); 228 229 PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC,Mat); 230 PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC,Mat); 231 PETSC_EXTERN PetscErrorCode PCGalerkinSetComputeSubmatrix(PC,PetscErrorCode (*)(PC,Mat,Mat,Mat*,void*),void*); 232 233 PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC,PetscInt,PetscInt,PetscReal*); 234 235 PETSC_EXTERN PetscErrorCode PCPythonSetType(PC,const char[]); 236 237 PETSC_EXTERN PetscErrorCode PCSetDM(PC,DM); 238 PETSC_EXTERN PetscErrorCode PCGetDM(PC,DM*); 239 240 PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC,void*); 241 PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC,void*); 242 243 PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC,PCPARMSGlobalType); 244 PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC,PCPARMSLocalType); 245 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC,PetscReal,PetscInt); 246 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC,PetscInt); 247 PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC,PetscBool); 248 PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC,PetscInt,PetscInt,PetscInt); 249 250 PETSC_EXTERN PetscErrorCode PCGAMGSetType( PC,PCGAMGType); 251 PETSC_EXTERN PetscErrorCode PCGAMGGetType( PC,PCGAMGType*); 252 PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC,PetscInt); 253 PETSC_EXTERN PetscErrorCode PCGAMGSetRepartition(PC,PetscBool); 254 PETSC_EXTERN PetscErrorCode PCGAMGASMSetUseAggs(PC,PetscBool); 255 PETSC_EXTERN PetscErrorCode PCGAMGSetUseParallelCoarseGridSolve(PC,PetscBool); 256 PETSC_EXTERN PetscErrorCode PCGAMGSetSolverType(PC,char[],PetscInt); 257 PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC,PetscReal[],PetscInt); 258 PETSC_EXTERN PetscErrorCode PCGAMGSetThresholdScale(PC,PetscReal); 259 PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC,PetscInt); 260 PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC,PetscInt); 261 PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC,PetscInt); 262 PETSC_EXTERN PetscErrorCode PCGAMGSetSymGraph(PC,PetscBool); 263 PETSC_EXTERN PetscErrorCode PCGAMGSetSquareGraph(PC,PetscInt); 264 PETSC_EXTERN PetscErrorCode PCGAMGSetReuseInterpolation(PC,PetscBool); 265 PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void); 266 PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void); 267 PETSC_EXTERN PetscErrorCode PCGAMGRegister(PCGAMGType,PetscErrorCode (*)(PC)); 268 269 PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC,PCGAMGClassicalType); 270 PETSC_EXTERN PetscErrorCode PCGAMGClassicalGetType(PC,PCGAMGClassicalType*); 271 272 PETSC_EXTERN PetscErrorCode PCBDDCSetDiscreteGradient(PC,Mat,PetscInt,PetscInt,PetscBool,PetscBool); 273 PETSC_EXTERN PetscErrorCode PCBDDCSetDivergenceMat(PC,Mat,PetscBool,IS); 274 PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisMat(PC,Mat,PetscBool); 275 PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesIS(PC,IS); 276 PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC,IS); 277 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC,PetscInt); 278 PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC,PetscInt); 279 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC,IS); 280 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC,IS); 281 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC,IS*); 282 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC,IS*); 283 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC,IS); 284 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC,IS); 285 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC,IS*); 286 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC,IS*); 287 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC,PetscInt,IS[]); 288 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC,PetscInt,IS[]); 289 PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode); 290 PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC,PetscBool,const char*,Mat*,PC*); 291 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat,Vec,Vec); 292 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat,Vec,Vec); 293 294 PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC,PetscBool); 295 PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC,PetscScalar); 296 PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC,Vec); 297 298 PETSC_EXTERN PetscInt PetscMGLevelId; 299 PETSC_EXTERN PetscErrorCode PCMGSetType(PC,PCMGType); 300 PETSC_EXTERN PetscErrorCode PCMGGetType(PC,PCMGType*); 301 PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC,PetscInt,MPI_Comm*); 302 PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC,PetscInt*); 303 304 PETSC_EXTERN PetscErrorCode PCMGSetDistinctSmoothUp(PC); 305 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmooth(PC,PetscInt); 306 PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC,PCMGCycleType); 307 PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC,PetscInt,PCMGCycleType); 308 PETSC_DEPRECATED("Use PCMGSetCycleTypeOnLevel") PETSC_STATIC_INLINE PetscErrorCode PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt t) {return PCMGSetCycleTypeOnLevel(pc,l,(PCMGCycleType)t);} 309 PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC,PetscInt); 310 PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC,PCMGGalerkinType); 311 PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC,PCMGGalerkinType*); 312 313 PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC,PetscInt,Vec); 314 PETSC_EXTERN PetscErrorCode PCMGSetX(PC,PetscInt,Vec); 315 PETSC_EXTERN PetscErrorCode PCMGSetR(PC,PetscInt,Vec); 316 317 PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC,PetscInt,Mat); 318 PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC,PetscInt,Mat*); 319 PETSC_EXTERN PetscErrorCode PCMGSetInjection(PC,PetscInt,Mat); 320 PETSC_EXTERN PetscErrorCode PCMGGetInjection(PC,PetscInt,Mat*); 321 PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC,PetscInt,Mat); 322 PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC,PetscInt,Mat*); 323 PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC,PetscInt,Vec); 324 PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC,PetscInt,Vec*); 325 PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC,PetscInt,PetscErrorCode (*)(Mat,Vec,Vec,Vec),Mat); 326 PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat,Vec,Vec,Vec); 327 328 PETSC_EXTERN PetscErrorCode PCTelescopeGetSubcommType(PC,PetscSubcommType*); 329 PETSC_EXTERN PetscErrorCode PCTelescopeSetSubcommType(PC,PetscSubcommType); 330 PETSC_EXTERN PetscErrorCode PCTelescopeGetReductionFactor(PC,PetscInt*); 331 PETSC_EXTERN PetscErrorCode PCTelescopeSetReductionFactor(PC,PetscInt); 332 PETSC_EXTERN PetscErrorCode PCTelescopeGetIgnoreDM(PC,PetscBool*); 333 PETSC_EXTERN PetscErrorCode PCTelescopeSetIgnoreDM(PC,PetscBool); 334 PETSC_EXTERN PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC,PetscBool*); 335 PETSC_EXTERN PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC,PetscBool); 336 PETSC_EXTERN PetscErrorCode PCTelescopeGetDM(PC,DM*); 337 338 #endif /* __PETSCPC_H */ 339