1 /* 2 Preconditioner module. 3 */ 4 #pragma once 5 6 #include <petscmat.h> 7 #include <petscdmtypes.h> 8 #include <petscpctypes.h> 9 10 /* MANSEC = KSP */ 11 /* SUBMANSEC = PC */ 12 13 PETSC_EXTERN PetscErrorCode PCInitializePackage(void); 14 PETSC_EXTERN PetscErrorCode PCFinalizePackage(void); 15 16 /* 17 PCList contains the list of preconditioners currently registered 18 These are added with PCRegister() 19 */ 20 PETSC_EXTERN PetscFunctionList PCList; 21 22 /* Logging support */ 23 PETSC_EXTERN PetscClassId PC_CLASSID; 24 25 /* Arrays of names for options in implementation PCs */ 26 PETSC_EXTERN const char *const *const PCSides; 27 PETSC_EXTERN const char *const PCJacobiTypes[]; 28 PETSC_EXTERN const char *const PCASMTypes[]; 29 PETSC_EXTERN const char *const PCGASMTypes[]; 30 PETSC_EXTERN const char *const PCCompositeTypes[]; 31 PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[]; 32 PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[]; 33 PETSC_EXTERN const char *const PCPARMSGlobalTypes[]; 34 PETSC_EXTERN const char *const PCPARMSLocalTypes[]; 35 PETSC_EXTERN const char *const PCMGTypes[]; 36 PETSC_EXTERN const char *const PCMGCycleTypes[]; 37 PETSC_EXTERN const char *const PCMGGalerkinTypes[]; 38 PETSC_EXTERN const char *const PCMGCoarseSpaceTypes[]; 39 PETSC_EXTERN const char *const PCExoticTypes[]; 40 PETSC_EXTERN const char *const PCPatchConstructTypes[]; 41 PETSC_EXTERN const char *const PCDeflationTypes[]; 42 PETSC_EXTERN const char *const *const PCFailedReasons; 43 44 PETSC_EXTERN PetscErrorCode PCCreate(MPI_Comm, PC *); 45 PETSC_EXTERN PetscErrorCode PCSetType(PC, PCType); 46 PETSC_EXTERN PetscErrorCode PCGetType(PC, PCType *); 47 PETSC_EXTERN PetscErrorCode PCSetUp(PC); 48 49 PETSC_EXTERN PetscErrorCode PCSetKSPNestLevel(PC, PetscInt); 50 PETSC_EXTERN PetscErrorCode PCGetKSPNestLevel(PC, PetscInt *); 51 52 PETSC_EXTERN PetscErrorCode PCSetFailedReason(PC, PCFailedReason); 53 PETSC_EXTERN PetscErrorCode PCGetFailedReason(PC, PCFailedReason *); 54 PETSC_DEPRECATED_FUNCTION(3, 11, 0, "PCGetFailedReason()", ) static inline PetscErrorCode PCGetSetUpFailedReason(PC pc, PCFailedReason *reason) 55 { 56 return PCGetFailedReason(pc, reason); 57 } 58 PETSC_EXTERN PetscErrorCode PCReduceFailedReason(PC); 59 60 PETSC_EXTERN PetscErrorCode PCSetUpOnBlocks(PC); 61 PETSC_EXTERN PetscErrorCode PCApply(PC, Vec, Vec); 62 PETSC_EXTERN PetscErrorCode PCMatApply(PC, Mat, Mat); 63 PETSC_EXTERN PetscErrorCode PCApplySymmetricLeft(PC, Vec, Vec); 64 PETSC_EXTERN PetscErrorCode PCApplySymmetricRight(PC, Vec, Vec); 65 PETSC_EXTERN PetscErrorCode PCApplyBAorAB(PC, PCSide, Vec, Vec, Vec); 66 PETSC_EXTERN PetscErrorCode PCApplyTranspose(PC, Vec, Vec); 67 PETSC_EXTERN PetscErrorCode PCApplyTransposeExists(PC, PetscBool *); 68 PETSC_EXTERN PetscErrorCode PCApplyBAorABTranspose(PC, PCSide, Vec, Vec, Vec); 69 PETSC_EXTERN PetscErrorCode PCSetReusePreconditioner(PC, PetscBool); 70 PETSC_EXTERN PetscErrorCode PCGetReusePreconditioner(PC, PetscBool *); 71 PETSC_EXTERN PetscErrorCode PCSetErrorIfFailure(PC, PetscBool); 72 73 #define PC_FILE_CLASSID 1211222 74 75 PETSC_EXTERN PetscErrorCode PCApplyRichardson(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *); 76 PETSC_EXTERN PetscErrorCode PCApplyRichardsonExists(PC, PetscBool *); 77 PETSC_EXTERN PetscErrorCode PCSetUseAmat(PC, PetscBool); 78 PETSC_EXTERN PetscErrorCode PCGetUseAmat(PC, PetscBool *); 79 80 PETSC_EXTERN PetscErrorCode PCRegister(const char[], PetscErrorCode (*)(PC)); 81 82 PETSC_EXTERN PetscErrorCode PCReset(PC); 83 PETSC_EXTERN PetscErrorCode PCDestroy(PC *); 84 PETSC_EXTERN PetscErrorCode PCSetFromOptions(PC); 85 86 PETSC_EXTERN PetscErrorCode PCFactorGetMatrix(PC, Mat *); 87 PETSC_EXTERN PetscErrorCode PCSetModifySubMatrices(PC, PetscErrorCode (*)(PC, PetscInt, const IS[], const IS[], Mat[], void *), void *); 88 PETSC_EXTERN PetscErrorCode PCModifySubMatrices(PC, PetscInt, const IS[], const IS[], Mat[], void *); 89 90 PETSC_EXTERN PetscErrorCode PCSetOperators(PC, Mat, Mat); 91 PETSC_EXTERN PetscErrorCode PCGetOperators(PC, Mat *, Mat *); 92 PETSC_EXTERN PetscErrorCode PCGetOperatorsSet(PC, PetscBool *, PetscBool *); 93 94 PETSC_EXTERN PetscErrorCode PCView(PC, PetscViewer); 95 PETSC_EXTERN PetscErrorCode PCLoad(PC, PetscViewer); 96 PETSC_EXTERN PetscErrorCode PCViewFromOptions(PC, PetscObject, const char[]); 97 98 PETSC_EXTERN PetscErrorCode PCSetOptionsPrefix(PC, const char[]); 99 PETSC_EXTERN PetscErrorCode PCAppendOptionsPrefix(PC, const char[]); 100 PETSC_EXTERN PetscErrorCode PCGetOptionsPrefix(PC, const char *[]); 101 102 PETSC_EXTERN PetscErrorCode PCComputeOperator(PC, MatType, Mat *); 103 PETSC_DEPRECATED_FUNCTION(3, 12, 0, "PCComputeOperator()", ) static inline PetscErrorCode PCComputeExplicitOperator(PC A, Mat *B) 104 { 105 return PCComputeOperator(A, PETSC_NULLPTR, B); 106 } 107 108 /* 109 These are used to provide extra scaling of preconditioned 110 operator for time-stepping schemes like in SUNDIALS 111 */ 112 PETSC_EXTERN PetscErrorCode PCGetDiagonalScale(PC, PetscBool *); 113 PETSC_EXTERN PetscErrorCode PCDiagonalScaleLeft(PC, Vec, Vec); 114 PETSC_EXTERN PetscErrorCode PCDiagonalScaleRight(PC, Vec, Vec); 115 PETSC_EXTERN PetscErrorCode PCSetDiagonalScale(PC, Vec); 116 117 PETSC_EXTERN PetscErrorCode PCSetDM(PC, DM); 118 PETSC_EXTERN PetscErrorCode PCGetDM(PC, DM *); 119 120 PETSC_EXTERN PetscErrorCode PCGetInterpolations(PC, PetscInt *, Mat *[]); 121 PETSC_EXTERN PetscErrorCode PCGetCoarseOperators(PC, PetscInt *, Mat *[]); 122 123 PETSC_EXTERN PetscErrorCode PCSetCoordinates(PC, PetscInt, PetscInt, PetscReal[]); 124 125 PETSC_EXTERN PetscErrorCode PCSetApplicationContext(PC, void *); 126 PETSC_EXTERN PetscErrorCode PCGetApplicationContext(PC, void *); 127 128 /* ------------- options specific to particular preconditioners --------- */ 129 130 PETSC_EXTERN PetscErrorCode PCJacobiSetType(PC, PCJacobiType); 131 PETSC_EXTERN PetscErrorCode PCJacobiGetType(PC, PCJacobiType *); 132 PETSC_EXTERN PetscErrorCode PCJacobiSetUseAbs(PC, PetscBool); 133 PETSC_EXTERN PetscErrorCode PCJacobiGetUseAbs(PC, PetscBool *); 134 PETSC_EXTERN PetscErrorCode PCJacobiSetFixDiagonal(PC, PetscBool); 135 PETSC_EXTERN PetscErrorCode PCJacobiGetFixDiagonal(PC, PetscBool *); 136 PETSC_EXTERN PetscErrorCode PCJacobiGetDiagonal(PC, Vec, Vec); 137 PETSC_EXTERN PetscErrorCode PCJacobiSetRowl1Scale(PC, PetscReal); 138 PETSC_EXTERN PetscErrorCode PCJacobiGetRowl1Scale(PC, PetscReal *); 139 PETSC_EXTERN PetscErrorCode PCSORSetSymmetric(PC, MatSORType); 140 PETSC_EXTERN PetscErrorCode PCSORGetSymmetric(PC, MatSORType *); 141 PETSC_EXTERN PetscErrorCode PCSORSetOmega(PC, PetscReal); 142 PETSC_EXTERN PetscErrorCode PCSORGetOmega(PC, PetscReal *); 143 PETSC_EXTERN PetscErrorCode PCSORSetIterations(PC, PetscInt, PetscInt); 144 PETSC_EXTERN PetscErrorCode PCSORGetIterations(PC, PetscInt *, PetscInt *); 145 146 PETSC_EXTERN PetscErrorCode PCEisenstatSetOmega(PC, PetscReal); 147 PETSC_EXTERN PetscErrorCode PCEisenstatGetOmega(PC, PetscReal *); 148 PETSC_EXTERN PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC, PetscBool); 149 PETSC_EXTERN PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC, PetscBool *); 150 151 PETSC_EXTERN PetscErrorCode PCBJacobiSetTotalBlocks(PC, PetscInt, const PetscInt[]); 152 PETSC_EXTERN PetscErrorCode PCBJacobiGetTotalBlocks(PC, PetscInt *, const PetscInt *[]); 153 PETSC_EXTERN PetscErrorCode PCBJacobiSetLocalBlocks(PC, PetscInt, const PetscInt[]); 154 PETSC_EXTERN PetscErrorCode PCBJacobiGetLocalBlocks(PC, PetscInt *, const PetscInt *[]); 155 156 PETSC_EXTERN PetscErrorCode PCShellSetApply(PC, PetscErrorCode (*)(PC, Vec, Vec)); 157 PETSC_EXTERN PetscErrorCode PCShellSetMatApply(PC, PetscErrorCode (*)(PC, Mat, Mat)); 158 PETSC_EXTERN PetscErrorCode PCShellSetApplySymmetricLeft(PC, PetscErrorCode (*)(PC, Vec, Vec)); 159 PETSC_EXTERN PetscErrorCode PCShellSetApplySymmetricRight(PC, PetscErrorCode (*)(PC, Vec, Vec)); 160 PETSC_EXTERN PetscErrorCode PCShellSetApplyBA(PC, PetscErrorCode (*)(PC, PCSide, Vec, Vec, Vec)); 161 PETSC_EXTERN PetscErrorCode PCShellSetApplyTranspose(PC, PetscErrorCode (*)(PC, Vec, Vec)); 162 PETSC_EXTERN PetscErrorCode PCShellSetSetUp(PC, PetscErrorCode (*)(PC)); 163 PETSC_EXTERN PetscErrorCode PCShellSetApplyRichardson(PC, PetscErrorCode (*)(PC, Vec, Vec, Vec, PetscReal, PetscReal, PetscReal, PetscInt, PetscBool, PetscInt *, PCRichardsonConvergedReason *)); 164 PETSC_EXTERN PetscErrorCode PCShellSetView(PC, PetscErrorCode (*)(PC, PetscViewer)); 165 PETSC_EXTERN PetscErrorCode PCShellSetDestroy(PC, PetscErrorCode (*)(PC)); 166 PETSC_EXTERN PetscErrorCode PCShellSetContext(PC, void *); 167 PETSC_EXTERN PetscErrorCode PCShellGetContext(PC, void *); 168 PETSC_EXTERN PetscErrorCode PCShellSetName(PC, const char[]); 169 PETSC_EXTERN PetscErrorCode PCShellGetName(PC, const char *[]); 170 171 PETSC_EXTERN PetscErrorCode PCFactorSetZeroPivot(PC, PetscReal); 172 173 PETSC_EXTERN PetscErrorCode PCFactorSetShiftType(PC, MatFactorShiftType); 174 PETSC_EXTERN PetscErrorCode PCFactorSetShiftAmount(PC, PetscReal); 175 176 PETSC_EXTERN PetscErrorCode PCFactorSetMatSolverType(PC, MatSolverType); 177 PETSC_EXTERN PetscErrorCode PCFactorGetMatSolverType(PC, MatSolverType *); 178 PETSC_EXTERN PetscErrorCode PCFactorSetUpMatSolverType(PC); 179 PETSC_DEPRECATED_FUNCTION(3, 9, 0, "PCFactorSetMatSolverType()", ) static inline PetscErrorCode PCFactorSetMatSolverPackage(PC pc, MatSolverType stype) 180 { 181 return PCFactorSetMatSolverType(pc, stype); 182 } 183 PETSC_DEPRECATED_FUNCTION(3, 9, 0, "PCFactorGetMatSolverType()", ) static inline PetscErrorCode PCFactorGetMatSolverPackage(PC pc, MatSolverType *stype) 184 { 185 return PCFactorGetMatSolverType(pc, stype); 186 } 187 PETSC_DEPRECATED_FUNCTION(3, 9, 0, "PCFactorSetUpMatSolverType()", ) static inline PetscErrorCode PCFactorSetUpMatSolverPackage(PC pc) 188 { 189 return PCFactorSetUpMatSolverType(pc); 190 } 191 192 PETSC_EXTERN PetscErrorCode PCFactorSetFill(PC, PetscReal); 193 PETSC_EXTERN PetscErrorCode PCFactorSetColumnPivot(PC, PetscReal); 194 PETSC_EXTERN PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC, PetscReal); 195 196 PETSC_EXTERN PetscErrorCode PCFactorSetMatOrderingType(PC, MatOrderingType); 197 PETSC_EXTERN PetscErrorCode PCFactorSetReuseOrdering(PC, PetscBool); 198 PETSC_EXTERN PetscErrorCode PCFactorSetReuseFill(PC, PetscBool); 199 PETSC_EXTERN PetscErrorCode PCFactorSetUseInPlace(PC, PetscBool); 200 PETSC_EXTERN PetscErrorCode PCFactorGetUseInPlace(PC, PetscBool *); 201 PETSC_EXTERN PetscErrorCode PCFactorSetAllowDiagonalFill(PC, PetscBool); 202 PETSC_EXTERN PetscErrorCode PCFactorGetAllowDiagonalFill(PC, PetscBool *); 203 PETSC_EXTERN PetscErrorCode PCFactorSetPivotInBlocks(PC, PetscBool); 204 205 PETSC_EXTERN PetscErrorCode PCFactorSetLevels(PC, PetscInt); 206 PETSC_EXTERN PetscErrorCode PCFactorGetLevels(PC, PetscInt *); 207 PETSC_EXTERN PetscErrorCode PCFactorSetDropTolerance(PC, PetscReal, PetscReal, PetscInt); 208 PETSC_EXTERN PetscErrorCode PCFactorGetZeroPivot(PC, PetscReal *); 209 PETSC_EXTERN PetscErrorCode PCFactorGetShiftAmount(PC, PetscReal *); 210 PETSC_EXTERN PetscErrorCode PCFactorGetShiftType(PC, MatFactorShiftType *); 211 212 PETSC_EXTERN PetscErrorCode PCASMSetLocalSubdomains(PC, PetscInt, IS[], IS[]); 213 PETSC_EXTERN PetscErrorCode PCASMSetTotalSubdomains(PC, PetscInt, IS[], IS[]); 214 PETSC_EXTERN PetscErrorCode PCASMSetOverlap(PC, PetscInt); 215 PETSC_EXTERN PetscErrorCode PCASMSetDMSubdomains(PC, PetscBool); 216 PETSC_EXTERN PetscErrorCode PCASMGetDMSubdomains(PC, PetscBool *); 217 PETSC_EXTERN PetscErrorCode PCASMSetSortIndices(PC, PetscBool); 218 219 PETSC_EXTERN PetscErrorCode PCASMSetType(PC, PCASMType); 220 PETSC_EXTERN PetscErrorCode PCASMGetType(PC, PCASMType *); 221 PETSC_EXTERN PetscErrorCode PCASMSetLocalType(PC, PCCompositeType); 222 PETSC_EXTERN PetscErrorCode PCASMGetLocalType(PC, PCCompositeType *); 223 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains(Mat, PetscInt, IS *[]); 224 PETSC_EXTERN PetscErrorCode PCASMDestroySubdomains(PetscInt, IS *[], IS *[]); 225 PETSC_EXTERN PetscErrorCode PCASMCreateSubdomains2D(PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, IS *[], IS *[]); 226 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubdomains(PC, PetscInt *, IS *[], IS *[]); 227 PETSC_EXTERN PetscErrorCode PCASMGetLocalSubmatrices(PC, PetscInt *, Mat *[]); 228 PETSC_EXTERN PetscErrorCode PCASMGetSubMatType(PC, MatType *); 229 PETSC_EXTERN PetscErrorCode PCASMSetSubMatType(PC, MatType); 230 231 PETSC_EXTERN PetscErrorCode PCGASMSetTotalSubdomains(PC, PetscInt); 232 PETSC_EXTERN PetscErrorCode PCGASMSetSubdomains(PC, PetscInt, IS[], IS[]); 233 PETSC_EXTERN PetscErrorCode PCGASMSetOverlap(PC, PetscInt); 234 PETSC_EXTERN PetscErrorCode PCGASMSetUseDMSubdomains(PC, PetscBool); 235 PETSC_EXTERN PetscErrorCode PCGASMGetUseDMSubdomains(PC, PetscBool *); 236 PETSC_EXTERN PetscErrorCode PCGASMSetSortIndices(PC, PetscBool); 237 238 PETSC_EXTERN PetscErrorCode PCGASMSetType(PC, PCGASMType); 239 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains(Mat, PetscInt, PetscInt *, IS *[]); 240 PETSC_EXTERN PetscErrorCode PCGASMDestroySubdomains(PetscInt, IS *[], IS *[]); 241 PETSC_EXTERN PetscErrorCode PCGASMCreateSubdomains2D(PC, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, IS *[], IS *[]); 242 PETSC_EXTERN PetscErrorCode PCGASMGetSubdomains(PC, PetscInt *, IS *[], IS *[]); 243 PETSC_EXTERN PetscErrorCode PCGASMGetSubmatrices(PC, PetscInt *, Mat *[]); 244 245 PETSC_EXTERN PetscErrorCode PCCompositeSetType(PC, PCCompositeType); 246 PETSC_EXTERN PetscErrorCode PCCompositeGetType(PC, PCCompositeType *); 247 PETSC_EXTERN PetscErrorCode PCCompositeAddPCType(PC, PCType); 248 PETSC_EXTERN PetscErrorCode PCCompositeAddPC(PC, PC); 249 PETSC_EXTERN PetscErrorCode PCCompositeGetNumberPC(PC, PetscInt *); 250 PETSC_EXTERN PetscErrorCode PCCompositeGetPC(PC, PetscInt, PC *); 251 PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlpha(PC, PetscScalar); 252 PETSC_EXTERN PetscErrorCode PCCompositeSpecialSetAlphaMat(PC, Mat); 253 254 PETSC_EXTERN PetscErrorCode PCRedundantSetNumber(PC, PetscInt); 255 PETSC_EXTERN PetscErrorCode PCRedundantSetScatter(PC, VecScatter, VecScatter); 256 PETSC_EXTERN PetscErrorCode PCRedundantGetOperators(PC, Mat *, Mat *); 257 258 PETSC_EXTERN PetscErrorCode PCSPAISetEpsilon(PC, PetscReal); 259 PETSC_EXTERN PetscErrorCode PCSPAISetNBSteps(PC, PetscInt); 260 PETSC_EXTERN PetscErrorCode PCSPAISetMax(PC, PetscInt); 261 PETSC_EXTERN PetscErrorCode PCSPAISetMaxNew(PC, PetscInt); 262 PETSC_EXTERN PetscErrorCode PCSPAISetBlockSize(PC, PetscInt); 263 PETSC_EXTERN PetscErrorCode PCSPAISetCacheSize(PC, PetscInt); 264 PETSC_EXTERN PetscErrorCode PCSPAISetVerbose(PC, PetscInt); 265 PETSC_EXTERN PetscErrorCode PCSPAISetSp(PC, PetscInt); 266 267 PETSC_EXTERN PetscErrorCode PCHYPRESetType(PC, const char[]); 268 PETSC_EXTERN PetscErrorCode PCHYPREGetType(PC, const char *[]); 269 PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteGradient(PC, Mat); 270 PETSC_EXTERN PetscErrorCode PCHYPRESetDiscreteCurl(PC, Mat); 271 PETSC_EXTERN PetscErrorCode PCHYPRESetInterpolations(PC, PetscInt, Mat, Mat[], Mat, Mat[]); 272 PETSC_EXTERN PetscErrorCode PCHYPRESetEdgeConstantVectors(PC, Vec, Vec, Vec); 273 PETSC_EXTERN PetscErrorCode PCHYPREAMSSetInteriorNodes(PC, Vec); 274 PETSC_EXTERN PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC, Mat); 275 PETSC_EXTERN PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC, Mat); 276 PETSC_EXTERN PetscErrorCode PCHYPREGetCFMarkers(PC pc, PetscInt *[], PetscBT *[]); 277 278 PETSC_EXTERN PetscErrorCode PCFieldSplitSetFields(PC, const char[], PetscInt, const PetscInt *, const PetscInt *); 279 PETSC_EXTERN PetscErrorCode PCFieldSplitSetType(PC, PCCompositeType); 280 PETSC_EXTERN PetscErrorCode PCFieldSplitGetType(PC, PCCompositeType *); 281 PETSC_EXTERN PetscErrorCode PCFieldSplitSetBlockSize(PC, PetscInt); 282 PETSC_EXTERN PetscErrorCode PCFieldSplitSetIS(PC, const char[], IS); 283 PETSC_EXTERN PetscErrorCode PCFieldSplitGetIS(PC, const char[], IS *); 284 PETSC_EXTERN PetscErrorCode PCFieldSplitGetISByIndex(PC, PetscInt, IS *); 285 PETSC_EXTERN PetscErrorCode PCFieldSplitRestrictIS(PC, IS); 286 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDMSplits(PC, PetscBool); 287 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDMSplits(PC, PetscBool *); 288 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDiagUseAmat(PC, PetscBool); 289 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDiagUseAmat(PC, PetscBool *); 290 PETSC_EXTERN PetscErrorCode PCFieldSplitSetOffDiagUseAmat(PC, PetscBool); 291 PETSC_EXTERN PetscErrorCode PCFieldSplitGetOffDiagUseAmat(PC, PetscBool *); 292 293 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "PCFieldSplitSetSchurPre()", ) PetscErrorCode PCFieldSplitSchurPrecondition(PC, PCFieldSplitSchurPreType, Mat); 294 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurPre(PC, PCFieldSplitSchurPreType, Mat); 295 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurPre(PC, PCFieldSplitSchurPreType *, Mat *); 296 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurFactType(PC, PCFieldSplitSchurFactType); 297 PETSC_EXTERN PetscErrorCode PCFieldSplitSetSchurScale(PC, PetscScalar); 298 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSchurBlocks(PC, Mat *, Mat *, Mat *, Mat *); 299 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetS(PC, Mat *); 300 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurRestoreS(PC, Mat *); 301 PETSC_EXTERN PetscErrorCode PCFieldSplitGetDetectSaddlePoint(PC, PetscBool *); 302 PETSC_EXTERN PetscErrorCode PCFieldSplitSetDetectSaddlePoint(PC, PetscBool); 303 PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBTol(PC, PetscReal); 304 PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBNu(PC, PetscReal); 305 PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBMaxit(PC, PetscInt); 306 PETSC_EXTERN PetscErrorCode PCFieldSplitSetGKBDelay(PC, PetscInt); 307 308 PETSC_EXTERN PetscErrorCode PCGalerkinSetRestriction(PC, Mat); 309 PETSC_EXTERN PetscErrorCode PCGalerkinSetInterpolation(PC, Mat); 310 PETSC_EXTERN PetscErrorCode PCGalerkinSetComputeSubmatrix(PC, PetscErrorCode (*)(PC, Mat, Mat, Mat *, void *), void *); 311 312 PETSC_EXTERN PetscErrorCode PCPythonSetType(PC, const char[]); 313 PETSC_EXTERN PetscErrorCode PCPythonGetType(PC, const char *[]); 314 315 PETSC_EXTERN PetscErrorCode PCPARMSSetGlobal(PC, PCPARMSGlobalType); 316 PETSC_EXTERN PetscErrorCode PCPARMSSetLocal(PC, PCPARMSLocalType); 317 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveTolerances(PC, PetscReal, PetscInt); 318 PETSC_EXTERN PetscErrorCode PCPARMSSetSolveRestart(PC, PetscInt); 319 PETSC_EXTERN PetscErrorCode PCPARMSSetNonsymPerm(PC, PetscBool); 320 PETSC_EXTERN PetscErrorCode PCPARMSSetFill(PC, PetscInt, PetscInt, PetscInt); 321 322 PETSC_EXTERN PetscErrorCode PCGAMGSetType(PC, PCGAMGType); 323 PETSC_EXTERN PetscErrorCode PCGAMGGetType(PC, PCGAMGType *); 324 PETSC_EXTERN PetscErrorCode PCGAMGSetProcEqLim(PC, PetscInt); 325 326 PETSC_EXTERN PetscErrorCode PCGAMGSetRepartition(PC, PetscBool); 327 PETSC_EXTERN PetscErrorCode PCGAMGSetUseSAEstEig(PC, PetscBool); 328 PETSC_EXTERN PetscErrorCode PCGAMGSetRecomputeEstEig(PC, PetscBool); 329 PETSC_EXTERN PetscErrorCode PCGAMGSetEigenvalues(PC, PetscReal, PetscReal); 330 PETSC_EXTERN PetscErrorCode PCGAMGASMSetUseAggs(PC, PetscBool); 331 PETSC_EXTERN PetscErrorCode PCGAMGSetParallelCoarseGridSolve(PC, PetscBool); 332 PETSC_EXTERN PetscErrorCode PCGAMGSetCpuPinCoarseGrids(PC, PetscBool); 333 PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseGridLayoutType(PC, PCGAMGLayoutType); 334 PETSC_EXTERN PetscErrorCode PCGAMGSetThreshold(PC, PetscReal[], PetscInt); 335 PETSC_EXTERN PetscErrorCode PCGAMGSetRankReductionFactors(PC, PetscInt[], PetscInt); 336 PETSC_EXTERN PetscErrorCode PCGAMGSetThresholdScale(PC, PetscReal); 337 PETSC_EXTERN PetscErrorCode PCGAMGSetCoarseEqLim(PC, PetscInt); 338 PETSC_EXTERN PetscErrorCode PCGAMGSetNlevels(PC, PetscInt); 339 PETSC_EXTERN PetscErrorCode PCGAMGSetNSmooths(PC, PetscInt); 340 PETSC_EXTERN PetscErrorCode PCGAMGSetAggressiveLevels(PC, PetscInt); 341 PETSC_EXTERN PetscErrorCode PCGAMGSetReuseInterpolation(PC, PetscBool); 342 PETSC_EXTERN PetscErrorCode PCGAMGFinalizePackage(void); 343 PETSC_EXTERN PetscErrorCode PCGAMGInitializePackage(void); 344 PETSC_EXTERN PetscErrorCode PCGAMGRegister(PCGAMGType, PetscErrorCode (*)(PC)); 345 PETSC_EXTERN PetscErrorCode PCGAMGCreateGraph(PC, Mat, Mat *); 346 PETSC_EXTERN PetscErrorCode PCGAMGSetAggressiveSquareGraph(PC, PetscBool); 347 PETSC_EXTERN PetscErrorCode PCGAMGMISkSetMinDegreeOrdering(PC, PetscBool); 348 PETSC_EXTERN PetscErrorCode PCGAMGMISkSetAggressive(PC, PetscInt); 349 PETSC_EXTERN PetscErrorCode PCGAMGASMSetHEM(PC, PetscInt); 350 PETSC_EXTERN PetscErrorCode PCGAMGSetLowMemoryFilter(PC, PetscBool); 351 PETSC_EXTERN PetscErrorCode PCGAMGSetGraphSymmetrize(PC, PetscBool); 352 PETSC_EXTERN PetscErrorCode PCGAMGSetInjectionIndex(PC, PetscInt, PetscInt[]); 353 354 PETSC_EXTERN PetscErrorCode PCGAMGClassicalSetType(PC, PCGAMGClassicalType); 355 PETSC_EXTERN PetscErrorCode PCGAMGClassicalGetType(PC, PCGAMGClassicalType *); 356 357 PETSC_EXTERN PetscErrorCode PCBDDCSetDiscreteGradient(PC, Mat, PetscInt, PetscInt, PetscBool, PetscBool); 358 PETSC_EXTERN PetscErrorCode PCBDDCSetDivergenceMat(PC, Mat, PetscBool, IS); 359 PETSC_EXTERN PetscErrorCode PCBDDCSetChangeOfBasisMat(PC, Mat, PetscBool); 360 PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesIS(PC, IS); 361 PETSC_EXTERN PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC, IS); 362 PETSC_EXTERN PetscErrorCode PCBDDCGetPrimalVerticesIS(PC, IS *); 363 PETSC_EXTERN PetscErrorCode PCBDDCGetPrimalVerticesLocalIS(PC, IS *); 364 PETSC_EXTERN PetscErrorCode PCBDDCSetCoarseningRatio(PC, PetscInt); 365 PETSC_EXTERN PetscErrorCode PCBDDCSetLevels(PC, PetscInt); 366 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundaries(PC, IS); 367 PETSC_EXTERN PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC, IS); 368 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundaries(PC, IS *); 369 PETSC_EXTERN PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC, IS *); 370 PETSC_EXTERN PetscErrorCode PCBDDCSetInterfaceExtType(PC, PCBDDCInterfaceExtType); 371 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundaries(PC, IS); 372 PETSC_EXTERN PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC, IS); 373 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundaries(PC, IS *); 374 PETSC_EXTERN PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC, IS *); 375 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplitting(PC, PetscInt, IS[]); 376 PETSC_EXTERN PetscErrorCode PCBDDCSetDofsSplittingLocal(PC, PetscInt, IS[]); 377 PETSC_EXTERN PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC, PetscInt, const PetscInt[], const PetscInt[], PetscCopyMode); 378 PETSC_EXTERN PetscErrorCode PCBDDCCreateFETIDPOperators(PC, PetscBool, const char *, Mat *, PC *); 379 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat, Vec, Vec); 380 PETSC_EXTERN PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat, Vec, Vec); 381 PETSC_EXTERN PetscErrorCode PCBDDCFinalizePackage(void); 382 PETSC_EXTERN PetscErrorCode PCBDDCInitializePackage(void); 383 384 PETSC_EXTERN PetscErrorCode PCISInitialize(PC); 385 PETSC_EXTERN PetscErrorCode PCISSetUp(PC, PetscBool, PetscBool); 386 PETSC_EXTERN PetscErrorCode PCISSetUseStiffnessScaling(PC, PetscBool); 387 PETSC_EXTERN PetscErrorCode PCISSetSubdomainScalingFactor(PC, PetscScalar); 388 PETSC_EXTERN PetscErrorCode PCISSetSubdomainDiagonalScaling(PC, Vec); 389 PETSC_EXTERN PetscErrorCode PCISScatterArrayNToVecB(PC, PetscScalar *, Vec, InsertMode, ScatterMode); 390 PETSC_EXTERN PetscErrorCode PCISApplySchur(PC, Vec, Vec, Vec, Vec, Vec); 391 PETSC_EXTERN PetscErrorCode PCISApplyInvSchur(PC, Vec, Vec, Vec, Vec); 392 PETSC_EXTERN PetscErrorCode PCISReset(PC); 393 394 PETSC_EXTERN PetscInt PetscMGLevelId; 395 PETSC_EXTERN PetscErrorCode PCMGSetType(PC, PCMGType); 396 PETSC_EXTERN PetscErrorCode PCMGGetType(PC, PCMGType *); 397 PETSC_EXTERN PetscErrorCode PCMGSetLevels(PC, PetscInt, MPI_Comm *); 398 PETSC_EXTERN PetscErrorCode PCMGGetLevels(PC, PetscInt *); 399 400 PETSC_EXTERN PetscErrorCode PCMGSetDistinctSmoothUp(PC); 401 PETSC_EXTERN PetscErrorCode PCMGSetNumberSmooth(PC, PetscInt); 402 PETSC_EXTERN PetscErrorCode PCMGSetCycleType(PC, PCMGCycleType); 403 PETSC_EXTERN PetscErrorCode PCMGSetCycleTypeOnLevel(PC, PetscInt, PCMGCycleType); 404 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "PCMGSetCycleTypeOnLevel()", ) static inline PetscErrorCode PCMGSetCyclesOnLevel(PC pc, PetscInt l, PetscInt t) 405 { 406 return PCMGSetCycleTypeOnLevel(pc, l, (PCMGCycleType)t); 407 } 408 PETSC_EXTERN PetscErrorCode PCMGMultiplicativeSetCycles(PC, PetscInt); 409 PETSC_EXTERN PetscErrorCode PCMGSetGalerkin(PC, PCMGGalerkinType); 410 PETSC_EXTERN PetscErrorCode PCMGGetGalerkin(PC, PCMGGalerkinType *); 411 PETSC_EXTERN PetscErrorCode PCMGSetAdaptCoarseSpaceType(PC, PCMGCoarseSpaceType); 412 PETSC_EXTERN PetscErrorCode PCMGGetAdaptCoarseSpaceType(PC, PCMGCoarseSpaceType *); 413 PETSC_EXTERN PetscErrorCode PCMGSetAdaptCR(PC, PetscBool); 414 PETSC_EXTERN PetscErrorCode PCMGGetAdaptCR(PC, PetscBool *); 415 /* MATT: Remove? */ 416 PETSC_EXTERN PetscErrorCode PCMGSetAdaptInterpolation(PC, PetscBool); 417 PETSC_EXTERN PetscErrorCode PCMGGetAdaptInterpolation(PC, PetscBool *); 418 419 PETSC_EXTERN PetscErrorCode PCMGSetRhs(PC, PetscInt, Vec); 420 PETSC_EXTERN PetscErrorCode PCMGSetX(PC, PetscInt, Vec); 421 PETSC_EXTERN PetscErrorCode PCMGSetR(PC, PetscInt, Vec); 422 423 PETSC_EXTERN PetscErrorCode PCMGSetRestriction(PC, PetscInt, Mat); 424 PETSC_EXTERN PetscErrorCode PCMGGetRestriction(PC, PetscInt, Mat *); 425 PETSC_EXTERN PetscErrorCode PCMGSetInjection(PC, PetscInt, Mat); 426 PETSC_EXTERN PetscErrorCode PCMGGetInjection(PC, PetscInt, Mat *); 427 PETSC_EXTERN PetscErrorCode PCMGSetInterpolation(PC, PetscInt, Mat); 428 PETSC_EXTERN PetscErrorCode PCMGSetOperators(PC, PetscInt, Mat, Mat); 429 PETSC_EXTERN PetscErrorCode PCMGGetInterpolation(PC, PetscInt, Mat *); 430 PETSC_EXTERN PetscErrorCode PCMGSetRScale(PC, PetscInt, Vec); 431 PETSC_EXTERN PetscErrorCode PCMGGetRScale(PC, PetscInt, Vec *); 432 PETSC_EXTERN PetscErrorCode PCMGSetResidual(PC, PetscInt, PetscErrorCode (*)(Mat, Vec, Vec, Vec), Mat); 433 PETSC_EXTERN PetscErrorCode PCMGSetResidualTranspose(PC, PetscInt, PetscErrorCode (*)(Mat, Vec, Vec, Vec), Mat); 434 PETSC_EXTERN PetscErrorCode PCMGResidualDefault(Mat, Vec, Vec, Vec); 435 PETSC_EXTERN PetscErrorCode PCMGResidualTransposeDefault(Mat, Vec, Vec, Vec); 436 PETSC_EXTERN PetscErrorCode PCMGMatResidualDefault(Mat, Mat, Mat, Mat); 437 PETSC_EXTERN PetscErrorCode PCMGMatResidualTransposeDefault(Mat, Mat, Mat, Mat); 438 PETSC_EXTERN PetscErrorCode PCMGGalerkinSetMatProductAlgorithm(PC, const char[]); 439 PETSC_EXTERN PetscErrorCode PCMGGalerkinGetMatProductAlgorithm(PC, const char *[]); 440 PETSC_EXTERN PetscErrorCode PCMGGetGridComplexity(PC, PetscReal *, PetscReal *); 441 442 PETSC_EXTERN PetscErrorCode PCHMGSetReuseInterpolation(PC, PetscBool); 443 PETSC_EXTERN PetscErrorCode PCHMGSetUseSubspaceCoarsening(PC, PetscBool); 444 PETSC_EXTERN PetscErrorCode PCHMGSetInnerPCType(PC, PCType); 445 PETSC_EXTERN PetscErrorCode PCHMGSetCoarseningComponent(PC, PetscInt); 446 PETSC_EXTERN PetscErrorCode PCHMGUseMatMAIJ(PC, PetscBool); 447 448 PETSC_EXTERN PetscErrorCode PCTelescopeGetSubcommType(PC, PetscSubcommType *); 449 PETSC_EXTERN PetscErrorCode PCTelescopeSetSubcommType(PC, PetscSubcommType); 450 PETSC_EXTERN PetscErrorCode PCTelescopeGetReductionFactor(PC, PetscInt *); 451 PETSC_EXTERN PetscErrorCode PCTelescopeSetReductionFactor(PC, PetscInt); 452 PETSC_EXTERN PetscErrorCode PCTelescopeGetIgnoreDM(PC, PetscBool *); 453 PETSC_EXTERN PetscErrorCode PCTelescopeSetIgnoreDM(PC, PetscBool); 454 PETSC_EXTERN PetscErrorCode PCTelescopeGetUseCoarseDM(PC, PetscBool *); 455 PETSC_EXTERN PetscErrorCode PCTelescopeSetUseCoarseDM(PC, PetscBool); 456 PETSC_EXTERN PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC, PetscBool *); 457 PETSC_EXTERN PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC, PetscBool); 458 PETSC_EXTERN PetscErrorCode PCTelescopeGetDM(PC, DM *); 459 460 PETSC_EXTERN PetscErrorCode PCPatchSetSaveOperators(PC, PetscBool); 461 PETSC_EXTERN PetscErrorCode PCPatchGetSaveOperators(PC, PetscBool *); 462 PETSC_EXTERN PetscErrorCode PCPatchSetPrecomputeElementTensors(PC, PetscBool); 463 PETSC_EXTERN PetscErrorCode PCPatchGetPrecomputeElementTensors(PC, PetscBool *); 464 PETSC_EXTERN PetscErrorCode PCPatchSetPartitionOfUnity(PC, PetscBool); 465 PETSC_EXTERN PetscErrorCode PCPatchGetPartitionOfUnity(PC, PetscBool *); 466 PETSC_EXTERN PetscErrorCode PCPatchSetSubMatType(PC, MatType); 467 PETSC_EXTERN PetscErrorCode PCPatchGetSubMatType(PC, MatType *); 468 PETSC_EXTERN PetscErrorCode PCPatchSetCellNumbering(PC, PetscSection); 469 PETSC_EXTERN PetscErrorCode PCPatchGetCellNumbering(PC, PetscSection *); 470 PETSC_EXTERN PetscErrorCode PCPatchSetConstructType(PC, PCPatchConstructType, PetscErrorCode (*)(PC, PetscInt *, IS **, IS *, void *), void *); 471 PETSC_EXTERN PetscErrorCode PCPatchGetConstructType(PC, PCPatchConstructType *, PetscErrorCode (**)(PC, PetscInt *, IS **, IS *, void *), void **); 472 PETSC_EXTERN PetscErrorCode PCPatchSetDiscretisationInfo(PC, PetscInt, DM *, PetscInt *, PetscInt *, const PetscInt **, const PetscInt *, PetscInt, const PetscInt *, PetscInt, const PetscInt *); 473 PETSC_EXTERN PetscErrorCode PCPatchSetComputeOperator(PC, PetscErrorCode (*)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *); 474 PETSC_EXTERN PetscErrorCode PCPatchSetComputeFunction(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx); 475 PETSC_EXTERN PetscErrorCode PCPatchSetComputeOperatorInteriorFacets(PC, PetscErrorCode (*)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *); 476 PETSC_EXTERN PetscErrorCode PCPatchSetComputeFunctionInteriorFacets(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx); 477 478 PETSC_EXTERN PetscErrorCode PCLMVMSetMatLMVM(PC, Mat); 479 PETSC_EXTERN PetscErrorCode PCLMVMGetMatLMVM(PC, Mat *); 480 PETSC_EXTERN PetscErrorCode PCLMVMSetIS(PC, IS); 481 PETSC_EXTERN PetscErrorCode PCLMVMClearIS(PC); 482 PETSC_EXTERN PetscErrorCode PCLMVMSetUpdateVec(PC, Vec); 483 484 PETSC_EXTERN PetscErrorCode PCExoticSetType(PC, PCExoticType); 485 486 PETSC_EXTERN PetscErrorCode PCDeflationSetInitOnly(PC, PetscBool); 487 PETSC_EXTERN PetscErrorCode PCDeflationSetLevels(PC, PetscInt); 488 PETSC_EXTERN PetscErrorCode PCDeflationSetReductionFactor(PC, PetscInt); 489 PETSC_EXTERN PetscErrorCode PCDeflationSetCorrectionFactor(PC, PetscScalar); 490 PETSC_EXTERN PetscErrorCode PCDeflationSetSpaceToCompute(PC, PCDeflationSpaceType, PetscInt); 491 PETSC_EXTERN PetscErrorCode PCDeflationSetSpace(PC, Mat, PetscBool); 492 PETSC_EXTERN PetscErrorCode PCDeflationSetProjectionNullSpaceMat(PC, Mat); 493 PETSC_EXTERN PetscErrorCode PCDeflationSetCoarseMat(PC, Mat); 494 PETSC_EXTERN PetscErrorCode PCDeflationGetPC(PC, PC *); 495 496 PETSC_EXTERN PetscErrorCode PCHPDDMSetAuxiliaryMat(PC, IS, Mat, PetscErrorCode (*)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void *), void *); 497 PETSC_EXTERN PetscErrorCode PCHPDDMSetRHSMat(PC, Mat); 498 PETSC_EXTERN PetscErrorCode PCHPDDMHasNeumannMat(PC, PetscBool); 499 PETSC_EXTERN PetscErrorCode PCHPDDMSetCoarseCorrectionType(PC, PCHPDDMCoarseCorrectionType); 500 PETSC_EXTERN PetscErrorCode PCHPDDMGetCoarseCorrectionType(PC, PCHPDDMCoarseCorrectionType *); 501 PETSC_EXTERN PetscErrorCode PCHPDDMSetSTShareSubKSP(PC, PetscBool); 502 PETSC_EXTERN PetscErrorCode PCHPDDMGetSTShareSubKSP(PC, PetscBool *); 503 PETSC_EXTERN PetscErrorCode PCHPDDMSetDeflationMat(PC, IS, Mat); 504 PETSC_EXTERN PetscErrorCode PCHPDDMFinalizePackage(void); 505 PETSC_EXTERN PetscErrorCode PCHPDDMInitializePackage(void); 506 PETSC_EXTERN PetscErrorCode PCHPDDMGetComplexities(PC, PetscReal *, PetscReal *); 507 508 PETSC_EXTERN PetscErrorCode PCAmgXGetResources(PC, void *); 509 510 PETSC_EXTERN PetscErrorCode PCMatSetApplyOperation(PC, MatOperation); 511 PETSC_EXTERN PetscErrorCode PCMatGetApplyOperation(PC, MatOperation *); 512