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(const char[]); 9 10 /*S 11 KSP - Abstract PETSc object that manages all Krylov methods 12 13 Level: beginner 14 15 Concepts: Krylov methods 16 17 .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP 18 S*/ 19 typedef struct _p_KSP* KSP; 20 21 /*J 22 KSPType - String with the name of a PETSc Krylov method or the creation function 23 with an optional dynamic library name, for example 24 http://www.mcs.anl.gov/petsc/lib.a:mykspcreate() 25 26 Level: beginner 27 28 .seealso: KSPSetType(), KSP 29 J*/ 30 typedef const char* KSPType; 31 #define KSPRICHARDSON "richardson" 32 #define KSPCHEBYSHEV "chebyshev" 33 #define KSPCG "cg" 34 #define KSPGROPPCG "groppcg" 35 #define KSPPIPECG "pipecg" 36 #define KSPCGNE "cgne" 37 #define KSPNASH "nash" 38 #define KSPSTCG "stcg" 39 #define KSPGLTR "gltr" 40 #define KSPGMRES "gmres" 41 #define KSPFGMRES "fgmres" 42 #define KSPLGMRES "lgmres" 43 #define KSPDGMRES "dgmres" 44 #define KSPPGMRES "pgmres" 45 #define KSPTCQMR "tcqmr" 46 #define KSPBCGS "bcgs" 47 #define KSPIBCGS "ibcgs" 48 #define KSPFBCGS "fbcgs" 49 #define KSPFBCGSR "fbcgsr" 50 #define KSPBCGSL "bcgsl" 51 #define KSPCGS "cgs" 52 #define KSPTFQMR "tfqmr" 53 #define KSPCR "cr" 54 #define KSPPIPECR "pipecr" 55 #define KSPLSQR "lsqr" 56 #define KSPPREONLY "preonly" 57 #define KSPQCG "qcg" 58 #define KSPBICG "bicg" 59 #define KSPMINRES "minres" 60 #define KSPSYMMLQ "symmlq" 61 #define KSPLCD "lcd" 62 #define KSPPYTHON "python" 63 #define KSPGCR "gcr" 64 #define KSPSPECEST "specest" 65 66 /* Logging support */ 67 PETSC_EXTERN PetscClassId KSP_CLASSID; 68 PETSC_EXTERN PetscClassId DMKSP_CLASSID; 69 70 PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *); 71 PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 72 PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 73 PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 74 PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 75 PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 76 PETSC_EXTERN PetscErrorCode KSPReset(KSP); 77 PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 78 79 PETSC_EXTERN PetscFList KSPList; 80 PETSC_EXTERN PetscBool KSPRegisterAllCalled; 81 PETSC_EXTERN PetscErrorCode KSPRegisterAll(const char[]); 82 PETSC_EXTERN PetscErrorCode KSPRegisterDestroy(void); 83 PETSC_EXTERN PetscErrorCode KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP)); 84 PETSC_EXTERN PetscErrorCode KSPMatRegisterAll(const char[]); 85 86 /*MC 87 KSPRegisterDynamic - Adds a method to the Krylov subspace solver package. 88 89 Synopsis: 90 PetscErrorCode KSPRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(KSP)) 91 92 Not Collective 93 94 Input Parameters: 95 + name_solver - name of a new user-defined solver 96 . path - path (either absolute or relative) the library containing this solver 97 . name_create - name of routine to create method context 98 - routine_create - routine to create method context 99 100 Notes: 101 KSPRegisterDynamic() may be called multiple times to add several user-defined solvers. 102 103 If dynamic libraries are used, then the fourth input argument (routine_create) 104 is ignored. 105 106 Sample usage: 107 .vb 108 KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 109 "MySolverCreate",MySolverCreate); 110 .ve 111 112 Then, your solver can be chosen with the procedural interface via 113 $ KSPSetType(ksp,"my_solver") 114 or at runtime via the option 115 $ -ksp_type my_solver 116 117 Level: advanced 118 119 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 120 and others of the form ${any_environmental_variable} occuring in pathname will be 121 replaced with appropriate values. 122 If your function is not being put into a shared library then use KSPRegister() instead 123 124 .keywords: KSP, register 125 126 .seealso: KSPRegisterAll(), KSPRegisterDestroy() 127 128 M*/ 129 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 130 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 131 #else 132 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 133 #endif 134 135 PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *); 136 PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 137 PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 138 PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 139 PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 140 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool ); 141 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *); 142 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool ); 143 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *); 144 PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool ); 145 PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *); 146 PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *); 147 PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 148 PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *); 149 PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool ); 150 PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *); 151 PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *); 152 PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 153 PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 154 PETSC_EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace); 155 PETSC_EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*); 156 PETSC_EXTERN PetscErrorCode KSPGetVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 157 158 PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 159 PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 160 161 PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 162 PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**)); 163 PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 164 PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **); 165 PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 166 PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 167 168 /* not sure where to put this */ 169 PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 170 PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 171 PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 172 PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 173 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 174 175 PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 176 177 PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 178 PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 179 180 PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 181 PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 182 PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 183 PETSC_EXTERN PetscErrorCode KSPChebyshevSetEstimateEigenvalues(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 184 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 185 PETSC_EXTERN PetscErrorCode KSPChebyshevSetNewMatrix(KSP); 186 PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 187 PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 188 PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 189 190 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 191 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 192 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 193 194 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 195 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 196 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 197 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 198 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 199 200 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 201 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 202 203 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 204 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 205 PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 206 207 /*E 208 KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 209 210 Level: advanced 211 212 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 213 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 214 215 E*/ 216 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 217 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 218 /*MC 219 KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 220 221 Level: advanced 222 223 Note: Possible unstable, but the fastest to compute 224 225 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 226 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 227 KSPGMRESModifiedGramSchmidtOrthogonalization() 228 M*/ 229 230 /*MC 231 KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 232 iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 233 poor orthogonality. 234 235 Level: advanced 236 237 Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 238 estimate the orthogonality but is more stable. 239 240 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 241 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 242 KSPGMRESModifiedGramSchmidtOrthogonalization() 243 M*/ 244 245 /*MC 246 KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 247 248 Level: advanced 249 250 Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 251 but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 252 253 You should only use this if you absolutely know that the iterative refinement is needed. 254 255 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 256 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 257 KSPGMRESModifiedGramSchmidtOrthogonalization() 258 M*/ 259 260 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 261 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 262 263 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 264 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 265 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 266 267 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 268 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 269 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 270 271 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 272 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 273 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 274 275 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 276 PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 277 278 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 279 PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 280 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 281 PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 282 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 283 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 284 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 285 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 286 PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 287 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 288 PETSC_EXTERN PetscErrorCode KSPMonitorAMS(KSP,PetscInt,PetscReal,void*); 289 PETSC_EXTERN PetscErrorCode KSPMonitorAMSCreate(KSP,const char*,void**); 290 PETSC_EXTERN PetscErrorCode KSPMonitorAMSDestroy(void**); 291 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 292 293 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 294 PETSC_EXTERN PetscErrorCode KSPDefaultBuildSolution(KSP,Vec,Vec*); 295 PETSC_EXTERN PetscErrorCode KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 296 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 297 298 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure); 299 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 300 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 301 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 302 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 303 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 304 PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 305 PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 306 307 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 308 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 309 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 310 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 311 312 PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 313 PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 314 315 #define KSP_FILE_CLASSID 1211223 316 317 PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 318 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 319 320 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 321 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 322 323 /*E 324 KSPNormType - Norm that is passed in the Krylov convergence 325 test routines. 326 327 Level: advanced 328 329 Each solver only supports a subset of these and some may support different ones 330 depending on left or right preconditioning, see KSPSetPCSide() 331 332 Notes: this must match finclude/petscksp.h 333 334 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 335 KSPSetConvergenceTest(), KSPSetPCSide() 336 E*/ 337 typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 338 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 339 PETSC_EXTERN const char *const*const KSPNormTypes; 340 341 /*MC 342 KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 343 possibly save some computation but means the convergence test cannot 344 be based on a norm of a residual etc. 345 346 Level: advanced 347 348 Note: Some Krylov methods need to compute a residual norm and then this option is ignored 349 350 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 351 M*/ 352 353 /*MC 354 KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual and pass that to the 355 convergence test routine. 356 357 Level: advanced 358 359 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 360 M*/ 361 362 /*MC 363 KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 364 convergence test routine. 365 366 Level: advanced 367 368 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 369 M*/ 370 371 /*MC 372 KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 373 convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS 374 375 Level: advanced 376 377 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 378 M*/ 379 380 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 381 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 382 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 383 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 384 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool ); 385 386 /*E 387 KSPConvergedReason - reason a Krylov method was said to 388 have converged or diverged 389 390 Level: beginner 391 392 Notes: See KSPGetConvergedReason() for explanation of each value 393 394 Developer notes: this must match finclude/petscksp.h 395 396 The string versions of these are KSPConvergedReasons; if you change 397 any of the values here also change them that array of names. 398 399 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 400 E*/ 401 typedef enum {/* converged */ 402 KSP_CONVERGED_RTOL_NORMAL = 1, 403 KSP_CONVERGED_ATOL_NORMAL = 9, 404 KSP_CONVERGED_RTOL = 2, 405 KSP_CONVERGED_ATOL = 3, 406 KSP_CONVERGED_ITS = 4, 407 KSP_CONVERGED_CG_NEG_CURVE = 5, 408 KSP_CONVERGED_CG_CONSTRAINED = 6, 409 KSP_CONVERGED_STEP_LENGTH = 7, 410 KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 411 /* diverged */ 412 KSP_DIVERGED_NULL = -2, 413 KSP_DIVERGED_ITS = -3, 414 KSP_DIVERGED_DTOL = -4, 415 KSP_DIVERGED_BREAKDOWN = -5, 416 KSP_DIVERGED_BREAKDOWN_BICG = -6, 417 KSP_DIVERGED_NONSYMMETRIC = -7, 418 KSP_DIVERGED_INDEFINITE_PC = -8, 419 KSP_DIVERGED_NAN = -9, 420 KSP_DIVERGED_INDEFINITE_MAT = -10, 421 422 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 423 PETSC_EXTERN const char *const*KSPConvergedReasons; 424 425 /*MC 426 KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 427 428 Level: beginner 429 430 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 431 for left preconditioning it is the 2-norm of the preconditioned residual, and the 432 2-norm of the residual for right preconditioning 433 434 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 435 436 M*/ 437 438 /*MC 439 KSP_CONVERGED_ATOL - norm(r) <= atol 440 441 Level: beginner 442 443 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 444 for left preconditioning it is the 2-norm of the preconditioned residual, and the 445 2-norm of the residual for right preconditioning 446 447 Level: beginner 448 449 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 450 451 M*/ 452 453 /*MC 454 KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 455 456 Level: beginner 457 458 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 459 for left preconditioning it is the 2-norm of the preconditioned residual, and the 460 2-norm of the residual for right preconditioning 461 462 Level: beginner 463 464 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 465 466 M*/ 467 468 /*MC 469 KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 470 reached 471 472 Level: beginner 473 474 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 475 476 M*/ 477 478 /*MC 479 KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 480 the preconditioner is applied. Also used when the KSPSkipConverged() convergence 481 test routine is set in KSP. 482 483 484 Level: beginner 485 486 487 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 488 489 M*/ 490 491 /*MC 492 KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 493 method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 494 preconditioner. 495 496 Level: beginner 497 498 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 499 500 M*/ 501 502 /*MC 503 KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 504 method could not continue to enlarge the Krylov space. 505 506 507 Level: beginner 508 509 510 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 511 512 M*/ 513 514 /*MC 515 KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 516 symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 517 518 Level: beginner 519 520 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 521 522 M*/ 523 524 /*MC 525 KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 526 positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 527 be positive definite 528 529 Level: beginner 530 531 Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 532 the PCICC preconditioner to generate a positive definite preconditioner 533 534 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 535 536 M*/ 537 538 /*MC 539 KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 540 while the KSPSolve() is still running. 541 542 Level: beginner 543 544 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 545 546 M*/ 547 548 PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 549 PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 550 PETSC_EXTERN PetscErrorCode KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 551 PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 552 PETSC_EXTERN PetscErrorCode KSPDefaultConvergedDestroy(void *); 553 PETSC_EXTERN PetscErrorCode KSPDefaultConvergedCreate(void **); 554 PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUIRNorm(KSP); 555 PETSC_EXTERN PetscErrorCode KSPDefaultConvergedSetUMIRNorm(KSP); 556 PETSC_EXTERN PetscErrorCode KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 557 PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 558 559 PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 560 561 /*E 562 KSPCGType - Determines what type of CG to use 563 564 Level: beginner 565 566 .seealso: KSPCGSetType() 567 E*/ 568 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 569 PETSC_EXTERN const char *const KSPCGTypes[]; 570 571 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 572 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 573 574 PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 575 PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 576 PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 577 578 PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 579 PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 580 PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 581 582 PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 583 PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 584 PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 585 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 586 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 587 588 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 589 590 PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 591 PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 592 593 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 594 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 595 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormDestroy(PetscDrawLG*); 596 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 597 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 598 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG*); 599 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 600 601 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 602 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 603 604 /* see src/ksp/ksp/interface/iguess.c */ 605 typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 606 607 PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 608 PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 609 PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 610 PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 611 PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 612 PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 613 614 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 615 PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 616 PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 617 618 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 619 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 620 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 621 PETSC_EXTERN PetscErrorCode MatSchurComplementSet(Mat,Mat,Mat,Mat,Mat,Mat); 622 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdate(Mat,Mat,Mat,Mat,Mat,Mat,MatStructure); 623 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubmatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 624 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatReuse,Mat *); 625 626 PETSC_EXTERN PetscErrorCode MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat); 627 628 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 629 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 630 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 631 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 632 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 633 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 634 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 635 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 636 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,MatStructure*,void*),void*); 637 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,MatStructure*,void*),void*); 638 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 639 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 640 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 641 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 642 643 #endif 644