1 /* 2 Defines the interface functions for the Krylov subspace accelerators. 3 */ 4 #ifndef __PETSCKSP_H 5 #define __PETSCKSP_H 6 #include <petscpc.h> 7 8 PETSC_EXTERN PetscErrorCode KSPInitializePackage(void); 9 10 /*S 11 KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the 12 linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators). 13 14 Level: beginner 15 16 Concepts: Krylov methods 17 18 Notes: When a direct solver is used but no Krylov solver is used the KSP object is still used by with a 19 KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver). 20 21 .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy() 22 S*/ 23 typedef struct _p_KSP* KSP; 24 25 /*J 26 KSPType - String with the name of a PETSc Krylov method. 27 28 Level: beginner 29 30 .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions() 31 J*/ 32 typedef const char* KSPType; 33 #define KSPRICHARDSON "richardson" 34 #define KSPCHEBYSHEV "chebyshev" 35 #define KSPCG "cg" 36 #define KSPGROPPCG "groppcg" 37 #define KSPPIPECG "pipecg" 38 #define KSPCGNE "cgne" 39 #define KSPNASH "nash" 40 #define KSPSTCG "stcg" 41 #define KSPGLTR "gltr" 42 #define KSPFCG "fcg" 43 #define KSPPIPEFCG "pipefcg" 44 #define KSPGMRES "gmres" 45 #define KSPFGMRES "fgmres" 46 #define KSPLGMRES "lgmres" 47 #define KSPDGMRES "dgmres" 48 #define KSPPGMRES "pgmres" 49 #define KSPTCQMR "tcqmr" 50 #define KSPBCGS "bcgs" 51 #define KSPIBCGS "ibcgs" 52 #define KSPFBCGS "fbcgs" 53 #define KSPFBCGSR "fbcgsr" 54 #define KSPBCGSL "bcgsl" 55 #define KSPCGS "cgs" 56 #define KSPTFQMR "tfqmr" 57 #define KSPCR "cr" 58 #define KSPPIPECR "pipecr" 59 #define KSPLSQR "lsqr" 60 #define KSPPREONLY "preonly" 61 #define KSPQCG "qcg" 62 #define KSPBICG "bicg" 63 #define KSPMINRES "minres" 64 #define KSPSYMMLQ "symmlq" 65 #define KSPLCD "lcd" 66 #define KSPPYTHON "python" 67 #define KSPGCR "gcr" 68 #define KSPPIPEGCR "pipegcr" 69 #define KSPTSIRM "tsirm" 70 #define KSPCGLS "cgls" 71 72 /* Logging support */ 73 PETSC_EXTERN PetscClassId KSP_CLASSID; 74 PETSC_EXTERN PetscClassId DMKSP_CLASSID; 75 76 PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *); 77 PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 78 PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *); 79 PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 80 PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 81 PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 82 PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 83 PETSC_EXTERN PetscErrorCode KSPReset(KSP); 84 PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 85 PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 86 PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool); 87 88 PETSC_EXTERN PetscFunctionList KSPList; 89 PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 90 91 PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 92 PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 93 PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 94 PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 95 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool ); 96 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *); 97 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool ); 98 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *); 99 PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool ); 100 PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *); 101 PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 102 PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool ); 103 PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *); 104 PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool ); 105 PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *); 106 PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *); 107 PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *); 108 PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 109 PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 110 PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*); 111 PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 112 PETSC_DEPRECATED("Use KSPCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);} 113 114 PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 115 PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 116 117 PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 118 PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 119 120 PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 121 PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**)); 122 PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 123 PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **); 124 PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 125 PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 126 127 PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 128 PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *); 129 PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 130 PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 131 132 PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 133 PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 134 PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 135 PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 136 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 137 PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 138 PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 139 PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 140 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 141 PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 142 143 PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 144 PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 145 146 PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 147 PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 148 PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 149 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 150 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseRandom(KSP,PetscBool); 151 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 152 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*); 153 PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 154 PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *); 155 PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 156 PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt *,Vec[],PetscReal[],PetscReal[]); 157 158 /*E 159 160 KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 161 162 KSP_FCD_TRUNCATION uses all (up to mmax) stored directions 163 KSP_FCD_TRUNCATION_RESTART uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 164 165 Level: intermediate 166 .seealso : KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType() 167 168 E*/ 169 typedef enum {KSP_FCD_TRUNCATION,KSP_FCD_TRUNCATION_RESTART} KSPFCDTruncationType; 170 PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 171 172 #define KSPFCDGetNumOldDirections(ctx,i,mi) ({ \ 173 if(ctx->trunctype == KSP_FCD_TRUNCATION_RESTART){ \ 174 mi = ((i-1) % ctx->mmax)+1; \ 175 if (mi==1 && i!=1) ++(ctx->n_search_space_resets); \ 176 } else if (ctx->trunctype == KSP_FCD_TRUNCATION) \ 177 mi = ctx->mmax; \ 178 else { \ 179 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unrecognized Truncation Strategy");CHKERRQ(ierr); \ 180 } \ 181 }) 182 183 PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 184 PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 185 PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 186 PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 187 PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType); 188 PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*); 189 190 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt); 191 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*); 192 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt); 193 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*); 194 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType); 195 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*); 196 197 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt); 198 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*); 199 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt); 200 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*); 201 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType); 202 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*); 203 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool); 204 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*); 205 206 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 207 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 208 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 209 210 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 211 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 212 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 213 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 214 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 215 216 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 217 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 218 219 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 220 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 221 PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 222 223 /*E 224 KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 225 226 Level: advanced 227 228 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 229 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 230 231 E*/ 232 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 233 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 234 /*MC 235 KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 236 237 Level: advanced 238 239 Note: Possible unstable, but the fastest to compute 240 241 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 242 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 243 KSPGMRESModifiedGramSchmidtOrthogonalization() 244 M*/ 245 246 /*MC 247 KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 248 iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 249 poor orthogonality. 250 251 Level: advanced 252 253 Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 254 estimate the orthogonality but is more stable. 255 256 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 257 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 258 KSPGMRESModifiedGramSchmidtOrthogonalization() 259 M*/ 260 261 /*MC 262 KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 263 264 Level: advanced 265 266 Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 267 but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 268 269 You should only use this if you absolutely know that the iterative refinement is needed. 270 271 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 272 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 273 KSPGMRESModifiedGramSchmidtOrthogonalization() 274 M*/ 275 276 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 277 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 278 279 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 280 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 281 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 282 283 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 284 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 285 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 286 287 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 288 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 289 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 290 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 291 292 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 293 PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 294 295 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 296 PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 297 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 298 PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 299 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 300 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 301 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 302 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 303 PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 304 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 305 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 306 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 307 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 308 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 309 310 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 311 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 312 313 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 314 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 315 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 316 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 317 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 318 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 319 PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 320 PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 321 322 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 323 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 324 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 325 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 326 327 PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 328 PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 329 PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 330 PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer); 331 PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP); 332 333 #define KSP_FILE_CLASSID 1211223 334 335 PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 336 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 337 338 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 339 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 340 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 341 342 /*E 343 KSPNormType - Norm that is passed in the Krylov convergence 344 test routines. 345 346 Level: advanced 347 348 Each solver only supports a subset of these and some may support different ones 349 depending on left or right preconditioning, see KSPSetPCSide() 350 351 Notes: this must match petsc/finclude/petscksp.h 352 353 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 354 KSPSetConvergenceTest(), KSPSetPCSide() 355 E*/ 356 typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 357 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 358 PETSC_EXTERN const char *const*const KSPNormTypes; 359 360 /*MC 361 KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 362 possibly save some computation but means the convergence test cannot 363 be based on a norm of a residual etc. 364 365 Level: advanced 366 367 Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 368 369 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 370 M*/ 371 372 /*MC 373 KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 374 convergence test routine. 375 376 Level: advanced 377 378 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 379 M*/ 380 381 /*MC 382 KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 383 convergence test routine. 384 385 Level: advanced 386 387 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 388 M*/ 389 390 /*MC 391 KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 392 convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR 393 394 Level: advanced 395 396 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 397 M*/ 398 399 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 400 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 401 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 402 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 403 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 404 405 /*E 406 KSPConvergedReason - reason a Krylov method was said to have converged or diverged 407 408 Level: beginner 409 410 Notes: See KSPGetConvergedReason() for explanation of each value 411 412 Developer notes: this must match petsc/finclude/petscksp.h 413 414 The string versions of these are KSPConvergedReasons; if you change 415 any of the values here also change them that array of names. 416 417 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 418 E*/ 419 typedef enum {/* converged */ 420 KSP_CONVERGED_RTOL_NORMAL = 1, 421 KSP_CONVERGED_ATOL_NORMAL = 9, 422 KSP_CONVERGED_RTOL = 2, 423 KSP_CONVERGED_ATOL = 3, 424 KSP_CONVERGED_ITS = 4, 425 KSP_CONVERGED_CG_NEG_CURVE = 5, 426 KSP_CONVERGED_CG_CONSTRAINED = 6, 427 KSP_CONVERGED_STEP_LENGTH = 7, 428 KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 429 /* diverged */ 430 KSP_DIVERGED_NULL = -2, 431 KSP_DIVERGED_ITS = -3, 432 KSP_DIVERGED_DTOL = -4, 433 KSP_DIVERGED_BREAKDOWN = -5, 434 KSP_DIVERGED_BREAKDOWN_BICG = -6, 435 KSP_DIVERGED_NONSYMMETRIC = -7, 436 KSP_DIVERGED_INDEFINITE_PC = -8, 437 KSP_DIVERGED_NANORINF = -9, 438 KSP_DIVERGED_INDEFINITE_MAT = -10, 439 KSP_DIVERGED_PCSETUP_FAILED = -11, 440 441 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 442 PETSC_EXTERN const char *const*KSPConvergedReasons; 443 444 /*MC 445 KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 446 447 Level: beginner 448 449 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 450 for left preconditioning it is the 2-norm of the preconditioned residual, and the 451 2-norm of the residual for right preconditioning 452 453 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 454 455 M*/ 456 457 /*MC 458 KSP_CONVERGED_ATOL - norm(r) <= atol 459 460 Level: beginner 461 462 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 463 for left preconditioning it is the 2-norm of the preconditioned residual, and the 464 2-norm of the residual for right preconditioning 465 466 Level: beginner 467 468 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 469 470 M*/ 471 472 /*MC 473 KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 474 475 Level: beginner 476 477 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 478 for left preconditioning it is the 2-norm of the preconditioned residual, and the 479 2-norm of the residual for right preconditioning 480 481 Level: beginner 482 483 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 484 485 M*/ 486 487 /*MC 488 KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 489 reached 490 491 Level: beginner 492 493 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 494 495 M*/ 496 497 /*MC 498 KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 499 the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 500 test routine is set in KSP. 501 502 Level: beginner 503 504 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 505 506 M*/ 507 508 /*MC 509 KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 510 method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 511 preconditioner. 512 513 Level: beginner 514 515 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 516 517 M*/ 518 519 /*MC 520 KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 521 method could not continue to enlarge the Krylov space. 522 523 Level: beginner 524 525 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 526 527 M*/ 528 529 /*MC 530 KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 531 symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 532 533 Level: beginner 534 535 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 536 537 M*/ 538 539 /*MC 540 KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 541 positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 542 be positive definite 543 544 Level: beginner 545 546 Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 547 the PCICC preconditioner to generate a positive definite preconditioner 548 549 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 550 551 M*/ 552 553 /*MC 554 KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 555 while the KSPSolve() is still running. 556 557 Level: beginner 558 559 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 560 561 M*/ 562 563 PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 564 PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 565 PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 566 PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 567 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 568 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 569 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 570 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 571 PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 572 PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 573 574 PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 575 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 576 PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 577 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 578 PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 579 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 580 PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 581 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 582 PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 583 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 584 PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 585 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 586 587 PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 588 589 /*E 590 KSPCGType - Determines what type of CG to use 591 592 Level: beginner 593 594 .seealso: KSPCGSetType() 595 E*/ 596 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 597 PETSC_EXTERN const char *const KSPCGTypes[]; 598 599 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 600 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 601 602 PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 603 PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 604 PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 605 606 PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 607 PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 608 PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 609 610 PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 611 PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 612 PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 613 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 614 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 615 616 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 617 618 PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 619 PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 620 621 #include <petscdrawtypes.h> 622 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 623 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 624 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 625 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 626 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 627 628 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 629 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 630 631 /* see src/ksp/ksp/interface/iguess.c */ 632 typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 633 634 PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 635 PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 636 PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 637 PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 638 PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 639 PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 640 641 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 642 PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 643 PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 644 645 /*E 646 MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 647 648 Level: intermediate 649 650 .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 651 E*/ 652 typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType; 653 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 654 655 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 656 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 657 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 658 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 659 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 660 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 661 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 662 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 663 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 664 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 665 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *); 666 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 667 668 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 669 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 670 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 671 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 672 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 673 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 674 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 675 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 676 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 677 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 678 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 679 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 680 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 681 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 682 683 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 684 PETSC_EXTERN PetscErrorCode DMPlexProjectField(DM, Vec, void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode, Vec); 685 #endif 686