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