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 PETSC_EXTERN_CXX_BEGIN 8 9 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitializePackage(const char[]); 10 11 /*S 12 KSP - Abstract PETSc object that manages all Krylov methods 13 14 Level: beginner 15 16 Concepts: Krylov methods 17 18 .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP 19 S*/ 20 typedef struct _p_KSP* KSP; 21 22 /*E 23 KSPType - String with the name of a PETSc Krylov method or the creation function 24 with an optional dynamic library name, for example 25 http://www.mcs.anl.gov/petsc/lib.a:mykspcreate() 26 27 Level: beginner 28 29 .seealso: KSPSetType(), KSP 30 E*/ 31 #define KSPRICHARDSON "richardson" 32 #define KSPCHEBYCHEV "chebychev" 33 #define KSPCG "cg" 34 #define KSPCGNE "cgne" 35 #define KSPGMRES "gmres" 36 #define KSPFGMRES "fgmres" 37 #define KSPLGMRES "lgmres" 38 #define KSPTCQMR "tcqmr" 39 #define KSPBCGS "bcgs" 40 #define KSPBCGSL "bcgsl" 41 #define KSPCGS "cgs" 42 #define KSPTFQMR "tfqmr" 43 #define KSPCR "cr" 44 #define KSPLSQR "lsqr" 45 #define KSPPREONLY "preonly" 46 #define KSPQCG "qcg" 47 #define KSPBICG "bicg" 48 #define KSPMINRES "minres" 49 #define KSPSYMMLQ "symmlq" 50 #define KSPType const char* 51 52 /* Logging support */ 53 extern PetscCookie PETSCKSP_DLLEXPORT KSP_COOKIE; 54 extern PetscEvent KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve; 55 56 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate(MPI_Comm,KSP *); 57 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetType(KSP,KSPType); 58 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUp(KSP); 59 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetUpOnBlocks(KSP); 60 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolve(KSP,Vec,Vec); 61 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSolveTranspose(KSP,Vec,Vec); 62 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDestroy(KSP); 63 64 extern PetscFList KSPList; 65 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterAll(const char[]); 66 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegisterDestroy(void); 67 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP)); 68 69 /*MC 70 KSPRegisterDynamic - Adds a method to the Krylov subspace solver package. 71 72 Synopsis: 73 PetscErrorCode KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP)) 74 75 Not Collective 76 77 Input Parameters: 78 + name_solver - name of a new user-defined solver 79 . path - path (either absolute or relative) the library containing this solver 80 . name_create - name of routine to create method context 81 - routine_create - routine to create method context 82 83 Notes: 84 KSPRegisterDynamic() may be called multiple times to add several user-defined solvers. 85 86 If dynamic libraries are used, then the fourth input argument (routine_create) 87 is ignored. 88 89 Sample usage: 90 .vb 91 KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 92 "MySolverCreate",MySolverCreate); 93 .ve 94 95 Then, your solver can be chosen with the procedural interface via 96 $ KSPSetType(ksp,"my_solver") 97 or at runtime via the option 98 $ -ksp_type my_solver 99 100 Level: advanced 101 102 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 103 and others of the form ${any_environmental_variable} occuring in pathname will be 104 replaced with appropriate values. 105 If your function is not being put into a shared library then use KSPRegister() instead 106 107 .keywords: KSP, register 108 109 .seealso: KSPRegisterAll(), KSPRegisterDestroy() 110 111 M*/ 112 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 113 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 114 #else 115 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 116 #endif 117 118 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetType(KSP,KSPType *); 119 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPreconditionerSide(KSP,PCSide); 120 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPreconditionerSide(KSP,PCSide*); 121 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 122 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 123 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessNonzero(KSP,PetscTruth); 124 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessNonzero(KSP,PetscTruth *); 125 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetInitialGuessKnoll(KSP,PetscTruth); 126 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetInitialGuessKnoll(KSP,PetscTruth*); 127 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeEigenvalues(KSP,PetscTruth*); 128 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeEigenvalues(KSP,PetscTruth); 129 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetComputeSingularValues(KSP,PetscTruth*); 130 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetComputeSingularValues(KSP,PetscTruth); 131 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetRhs(KSP,Vec *); 132 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetSolution(KSP,Vec *); 133 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualNorm(KSP,PetscReal*); 134 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetIterationNumber(KSP,PetscInt*); 135 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNullSpace(KSP,MatNullSpace); 136 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetNullSpace(KSP,MatNullSpace*); 137 138 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetPC(KSP,PC); 139 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetPC(KSP,PC*); 140 141 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetMonitor(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void*)); 142 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPClearMonitor(KSP); 143 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetMonitorContext(KSP,void **); 144 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 145 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscTruth); 146 147 /* not sure where to put this */ 148 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCKSPGetKSP(PC,KSP*); 149 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 150 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 151 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 152 153 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildSolution(KSP,Vec,Vec *); 154 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPBuildResidual(KSP,Vec,Vec,Vec *); 155 156 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPRichardsonSetScale(KSP,PetscReal); 157 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal); 158 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 159 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvalues(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt *); 160 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal*,PetscReal*); 161 162 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart(KSP, PetscInt); 163 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetHapTol(KSP,PetscReal); 164 165 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors(KSP); 166 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 167 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 168 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 169 170 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetAugDim(KSP,PetscInt); 171 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMRESSetConstant(KSP); 172 173 /*E 174 KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 175 176 Level: advanced 177 178 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 179 KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 180 181 E*/ 182 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 183 extern const char *KSPGMRESCGSRefinementTypes[]; 184 /*M 185 KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 186 187 Level: advanced 188 189 Note: Possible unstable, but the fastest to compute 190 191 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 192 KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 193 KSPGMRESModifiedGramSchmidtOrthogonalization() 194 M*/ 195 196 /*M 197 KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 198 iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 199 poor orthogonality. 200 201 Level: advanced 202 203 Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 204 estimate the orthogonality but is more stable. 205 206 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 207 KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 208 KSPGMRESModifiedGramSchmidtOrthogonalization() 209 M*/ 210 211 /*M 212 KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 213 214 Level: advanced 215 216 Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 217 but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 218 219 You should only use this if you absolutely know that the iterative refinement is needed. 220 221 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 222 KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 223 KSPGMRESModifiedGramSchmidtOrthogonalization() 224 M*/ 225 226 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 227 228 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 229 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 230 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 231 232 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGSetTrustRegionRadius(KSP,PetscReal); 233 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetQuadratic(KSP,PetscReal*); 234 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPQCGGetTrialStepNorm(KSP,PetscReal*); 235 236 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetFromOptions(KSP); 237 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 238 239 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSingularValueMonitor(KSP,PetscInt,PetscReal,void *); 240 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultMonitor(KSP,PetscInt,PetscReal,void *); 241 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPTrueMonitor(KSP,PetscInt,PetscReal,void *); 242 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultSMonitor(KSP,PetscInt,PetscReal,void *); 243 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPVecViewMonitor(KSP,PetscInt,PetscReal,void *); 244 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESKrylovMonitor(KSP,PetscInt,PetscReal,void *); 245 246 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPUnwindPreconditioner(KSP,Vec,Vec); 247 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildSolution(KSP,Vec,Vec*); 248 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 249 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 250 251 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOperators(KSP,Mat,Mat,MatStructure); 252 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOperators(KSP,Mat*,Mat*,MatStructure*); 253 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetOptionsPrefix(KSP,const char[]); 254 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPAppendOptionsPrefix(KSP,const char[]); 255 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetOptionsPrefix(KSP,const char*[]); 256 257 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScale(KSP,PetscTruth); 258 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScale(KSP,PetscTruth*); 259 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetDiagonalScaleFix(KSP,PetscTruth); 260 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetDiagonalScaleFix(KSP,PetscTruth*); 261 262 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPView(KSP,PetscViewer); 263 264 /*E 265 KSPNormType - Norm that is passed in the Krylov convergence 266 test routines. 267 268 Level: advanced 269 270 Notes: this must match finclude/petscksp.h 271 272 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 273 KSPSetConvergenceTest() 274 E*/ 275 typedef enum {KSP_NO_NORM = 0, 276 KSP_PRECONDITIONED_NORM = 1, 277 KSP_UNPRECONDITIONED_NORM = 2, 278 KSP_NATURAL_NORM = 3} KSPNormType; 279 extern const char *KSPNormTypes[]; 280 /*M 281 KSP_NO_NORM - Do not compute a norm during the Krylov process. This will 282 possibly save some computation but means the convergence test cannot 283 be based on a norm of a residual etc. 284 285 Level: advanced 286 287 Note: Some Krylov methods need to compute a residual norm and then this option is ignored 288 289 .seealso: KSPNormType, KSPSetNormType(), KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM 290 M*/ 291 292 /*M 293 KSP_PRECONDITIONED_NORM - Compute the norm of the preconditioned residual and pass that to the 294 convergence test routine. 295 296 Level: advanced 297 298 .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest() 299 M*/ 300 301 /*M 302 KSP_UNPRECONDITIONED_NORM - Compute the norm of the true residual (b - A*x) and pass that to the 303 convergence test routine. 304 305 Level: advanced 306 307 .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest() 308 M*/ 309 310 /*M 311 KSP_NATURAL_NORM - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 312 convergence test routine. 313 314 Level: advanced 315 316 .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSPSetConvergenceTest() 317 M*/ 318 319 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetNormType(KSP,KSPNormType); 320 321 /*E 322 KSPConvergedReason - reason a Krylov method was said to 323 have converged or diverged 324 325 Level: beginner 326 327 Notes: this must match finclude/petscksp.h 328 329 Developer note: The string versions of these are in 330 src/ksp/ksp/interface/itfunc.c called convergedreasons. 331 If these enums are changed you must change those. 332 333 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 334 E*/ 335 typedef enum {/* converged */ 336 KSP_CONVERGED_RTOL = 2, 337 KSP_CONVERGED_ATOL = 3, 338 KSP_CONVERGED_ITS = 4, 339 KSP_CONVERGED_QCG_NEG_CURVE = 5, 340 KSP_CONVERGED_QCG_CONSTRAINED = 6, 341 KSP_CONVERGED_STEP_LENGTH = 7, 342 /* diverged */ 343 KSP_DIVERGED_NULL = -2, 344 KSP_DIVERGED_ITS = -3, 345 KSP_DIVERGED_DTOL = -4, 346 KSP_DIVERGED_BREAKDOWN = -5, 347 KSP_DIVERGED_BREAKDOWN_BICG = -6, 348 KSP_DIVERGED_NONSYMMETRIC = -7, 349 KSP_DIVERGED_INDEFINITE_PC = -8, 350 KSP_DIVERGED_NAN = -9, 351 KSP_DIVERGED_INDEFINITE_MAT = -10, 352 353 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 354 extern const char **KSPConvergedReasons; 355 356 /*MC 357 KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 358 359 Level: beginner 360 361 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 362 for left preconditioning it is the 2-norm of the preconditioned residual, and the 363 2-norm of the residual for right preconditioning 364 365 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 366 367 M*/ 368 369 /*MC 370 KSP_CONVERGED_ATOL - norm(r) <= atol 371 372 Level: beginner 373 374 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 375 for left preconditioning it is the 2-norm of the preconditioned residual, and the 376 2-norm of the residual for right preconditioning 377 378 Level: beginner 379 380 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 381 382 M*/ 383 384 /*MC 385 KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 386 387 Level: beginner 388 389 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 390 for left preconditioning it is the 2-norm of the preconditioned residual, and the 391 2-norm of the residual for right preconditioning 392 393 Level: beginner 394 395 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 396 397 M*/ 398 399 /*MC 400 KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 401 reached 402 403 Level: beginner 404 405 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 406 407 M*/ 408 409 /*MC 410 KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of the 411 preconditioner is applied. 412 413 414 Level: beginner 415 416 417 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 418 419 M*/ 420 421 /*MC 422 KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 423 method could not continue to enlarge the Krylov space. 424 425 Level: beginner 426 427 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 428 429 M*/ 430 431 /*MC 432 KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 433 method could not continue to enlarge the Krylov space. 434 435 436 Level: beginner 437 438 439 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 440 441 M*/ 442 443 /*MC 444 KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 445 symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 446 447 Level: beginner 448 449 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 450 451 M*/ 452 453 /*MC 454 KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 455 positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 456 be positive definite 457 458 Level: beginner 459 460 Notes: This can happen with the PCICC preconditioner, use -pc_icc_shift to force 461 the PCICC preconditioner to generate a positive definite preconditioner 462 463 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 464 465 M*/ 466 467 /*MC 468 KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 469 while the KSPSolve() is still running. 470 471 Level: beginner 472 473 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 474 475 M*/ 476 477 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *); 478 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergenceContext(KSP,void **); 479 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPDefaultConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 480 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPSkipConverged(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 481 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGetConvergedReason(KSP,KSPConvergedReason *); 482 483 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPComputeExplicitOperator(KSP,Mat *); 484 485 /*E 486 KSPCGType - Determines what type of CG to use 487 488 Level: beginner 489 490 .seealso: KSPCGSetType() 491 E*/ 492 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 493 extern const char *KSPCGTypes[]; 494 495 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPCGSetType(KSP,KSPCGType); 496 497 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPreSolve(PC,KSP); 498 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCPostSolve(PC,KSP); 499 500 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 501 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitor(KSP,PetscInt,PetscReal,void*); 502 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGMonitorDestroy(PetscDrawLG); 503 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 504 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitor(KSP,PetscInt,PetscReal,void*); 505 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPLGTrueMonitorDestroy(PetscDrawLG); 506 507 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostPreSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec)); 508 EXTERN PetscErrorCode PETSCKSP_DLLEXPORT PCShellSetPostSolve(PC,PetscErrorCode (*)(void*,KSP,Vec,Vec)); 509 510 PETSC_EXTERN_CXX_END 511 #endif 512