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