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