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 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 KSPTCQMR "tcqmr" 37 #define KSPBCGS "bcgs" 38 #define KSPBCGSL "bcgsl" 39 #define KSPCGS "cgs" 40 #define KSPTFQMR "tfqmr" 41 #define KSPCR "cr" 42 #define KSPLSQR "lsqr" 43 #define KSPPREONLY "preonly" 44 #define KSPQCG "qcg" 45 #define KSPBICG "bicg" 46 #define KSPFGMRES "fgmres" 47 #define KSPMINRES "minres" 48 #define KSPSYMMLQ "symmlq" 49 #define KSPLGMRES "lgmres" 50 #define KSPType char* 51 52 /* Logging support */ 53 extern PetscCookie KSP_COOKIE; 54 extern PetscEvent KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve; 55 56 EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *); 57 EXTERN PetscErrorCode KSPSetType(KSP,const KSPType); 58 EXTERN PetscErrorCode KSPSetUp(KSP); 59 EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 60 EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 61 EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 62 EXTERN PetscErrorCode KSPDestroy(KSP); 63 64 extern PetscFList KSPList; 65 EXTERN PetscErrorCode KSPRegisterAll(const char[]); 66 EXTERN PetscErrorCode KSPRegisterDestroy(void); 67 68 EXTERN PetscErrorCode KSPRegister(const char[],const char[],const char[],PetscErrorCode (*)(KSP)); 69 70 /*MC 71 KSPRegisterDynamic - Adds a method to the Krylov subspace solver package. 72 73 Synopsis: 74 int KSPRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(KSP)) 75 76 Not Collective 77 78 Input Parameters: 79 + name_solver - name of a new user-defined solver 80 . path - path (either absolute or relative) the library containing this solver 81 . name_create - name of routine to create method context 82 - routine_create - routine to create method context 83 84 Notes: 85 KSPRegisterDynamic() may be called multiple times to add several user-defined solvers. 86 87 If dynamic libraries are used, then the fourth input argument (routine_create) 88 is ignored. 89 90 Sample usage: 91 .vb 92 KSPRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 93 "MySolverCreate",MySolverCreate); 94 .ve 95 96 Then, your solver can be chosen with the procedural interface via 97 $ KSPSetType(ksp,"my_solver") 98 or at runtime via the option 99 $ -ksp_type my_solver 100 101 Level: advanced 102 103 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT}, 104 and others of the form ${any_environmental_variable} occuring in pathname will be 105 replaced with appropriate values. 106 If your function is not being put into a shared library then use KSPRegister() instead 107 108 .keywords: KSP, register 109 110 .seealso: KSPRegisterAll(), KSPRegisterDestroy() 111 112 M*/ 113 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 114 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,0) 115 #else 116 #define KSPRegisterDynamic(a,b,c,d) KSPRegister(a,b,c,d) 117 #endif 118 119 EXTERN PetscErrorCode KSPGetType(KSP,KSPType *); 120 EXTERN PetscErrorCode KSPSetPreconditionerSide(KSP,PCSide); 121 EXTERN PetscErrorCode KSPGetPreconditionerSide(KSP,PCSide*); 122 EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,int*); 123 EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,int); 124 EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscTruth); 125 EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscTruth *); 126 EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscTruth); 127 EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscTruth*); 128 EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscTruth); 129 EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscTruth); 130 EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *); 131 EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *); 132 EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 133 EXTERN PetscErrorCode KSPGetIterationNumber(KSP,int*); 134 EXTERN PetscErrorCode KSPSetNullSpace(KSP,MatNullSpace); 135 EXTERN PetscErrorCode KSPGetNullSpace(KSP,MatNullSpace*); 136 137 EXTERN PetscErrorCode KSPSetPC(KSP,PC); 138 EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 139 140 EXTERN PetscErrorCode KSPSetMonitor(KSP,PetscErrorCode (*)(KSP,int,PetscReal,void*),void *,PetscErrorCode (*)(void*)); 141 EXTERN PetscErrorCode KSPClearMonitor(KSP); 142 EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **); 143 EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],int *); 144 EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],int,PetscTruth); 145 146 /* not sure where to put this */ 147 EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 148 EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,int*,int*,KSP*[]); 149 EXTERN PetscErrorCode PCASMGetSubKSP(PC,int*,int*,KSP*[]); 150 151 EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 152 EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 153 154 EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 155 EXTERN PetscErrorCode KSPChebychevSetEigenvalues(KSP,PetscReal,PetscReal); 156 EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 157 EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,int,PetscReal*,PetscReal*,int *); 158 EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,int,PetscReal*,PetscReal*); 159 160 EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, int); 161 EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 162 163 EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 164 EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,int)); 165 EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,int); 166 EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,int); 167 168 EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,int); 169 EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 170 171 /*E 172 KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 173 174 Level: advanced 175 176 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 177 KSPGMRESSetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 178 179 E*/ 180 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 181 182 /*M 183 KSP_GMRES_CGS_REFINE_NEVER - Just do the classical (unmodified) Gram-Schmidt process 184 185 Level: advanced 186 187 Note: Possible unstable, but the fastest to compute 188 189 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 190 KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 191 KSPGMRESModifiedGramSchmidtOrthogonalization() 192 M*/ 193 194 /*M 195 KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 196 iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 197 poor orthogonality. 198 199 Level: advanced 200 201 Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 202 estimate the orthogonality but is more stable. 203 204 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 205 KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 206 KSPGMRESModifiedGramSchmidtOrthogonalization() 207 M*/ 208 209 /*M 210 KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 211 212 Level: advanced 213 214 Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 215 but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 216 217 You should only use this if you absolutely know that the iterative refinement is needed. 218 219 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), 220 KSPGMRESSetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 221 KSPGMRESModifiedGramSchmidtOrthogonalization() 222 M*/ 223 224 EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 225 226 EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,int,int,PetscReal,void*); 227 EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,int,int,PetscReal,void*); 228 EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,int,int,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 229 230 EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 231 EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 232 EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 233 234 EXTERN PetscErrorCode KSPSetFromOptions(KSP); 235 EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 236 237 EXTERN PetscErrorCode KSPSingularValueMonitor(KSP,int,PetscReal,void *); 238 EXTERN PetscErrorCode KSPDefaultMonitor(KSP,int,PetscReal,void *); 239 EXTERN PetscErrorCode KSPTrueMonitor(KSP,int,PetscReal,void *); 240 EXTERN PetscErrorCode KSPDefaultSMonitor(KSP,int,PetscReal,void *); 241 EXTERN PetscErrorCode KSPVecViewMonitor(KSP,int,PetscReal,void *); 242 EXTERN PetscErrorCode KSPGMRESKrylovMonitor(KSP,int,PetscReal,void *); 243 244 EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 245 EXTERN PetscErrorCode KSPDefaultBuildSolution(KSP,Vec,Vec*); 246 EXTERN PetscErrorCode KSPDefaultBuildResidual(KSP,Vec,Vec,Vec *); 247 248 EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat,MatStructure); 249 EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 250 EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 251 EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,char*[]); 252 253 EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscTruth); 254 EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscTruth*); 255 EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscTruth); 256 EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscTruth*); 257 258 EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 259 260 /*E 261 KSPNormType - Norm that is passed in the Krylov convergence 262 test routines. 263 264 Level: advanced 265 266 Notes: this must match finclude/petscksp.h 267 268 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 269 KSPSetConvergenceTest() 270 E*/ 271 typedef enum {KSP_NO_NORM = 0, 272 KSP_PRECONDITIONED_NORM = 1, 273 KSP_UNPRECONDITIONED_NORM = 2, 274 KSP_NATURAL_NORM = 3} KSPNormType; 275 276 /*M 277 KSP_NO_NORM - Do not compute a norm during the Krylov process. This will 278 possibly save some computation but means the convergence test cannot 279 be based on a norm of a residual etc. 280 281 Level: advanced 282 283 Note: Some Krylov methods need to compute a residual norm and then this is ignored 284 285 .seealso: KSPNormType, KSPSetNormType(), KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM 286 M*/ 287 288 /*M 289 KSP_PRECONDITIONED_NORM - Compute the norm of the preconditioned residual and pass that to the 290 convergence test routine. 291 292 Level: advanced 293 294 .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_UNPRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest() 295 M*/ 296 297 /*M 298 KSP_UNPRECONDITIONED_NORM - Compute the norm of the true residual (b - A*x) and pass that to the 299 convergence test routine. 300 301 Level: advanced 302 303 .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_NATURAL_NORM, KSPSetConvergenceTest() 304 M*/ 305 306 /*M 307 KSP_NATURAL_NORM - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 308 convergence test routine. 309 310 Level: advanced 311 312 .seealso: KSPNormType, KSPSetNormType(), KSP_NO_NORM, KSP_PRECONDITIONED_NORM, KSP_UNPRECONDITIONED_NORM, KSPSetConvergenceTest() 313 M*/ 314 315 EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 316 317 /*E 318 KSPConvergedReason - reason a Krylov method was said to 319 have converged or diverged 320 321 Level: beginner 322 323 Notes: this must match finclude/petscksp.h 324 325 Developer note: The string versions of these are in 326 src/ksp/ksp/interface/itfunc.c called convergedreasons. 327 If these enums are changed you much change those. 328 329 .seealso: KSPSolve(), KSPGetConvergedReason() 330 E*/ 331 typedef enum {/* converged */ 332 KSP_CONVERGED_RTOL = 2, 333 KSP_CONVERGED_ATOL = 3, 334 KSP_CONVERGED_ITS = 4, 335 KSP_CONVERGED_QCG_NEG_CURVE = 5, 336 KSP_CONVERGED_QCG_CONSTRAINED = 6, 337 KSP_CONVERGED_STEP_LENGTH = 7, 338 /* diverged */ 339 KSP_DIVERGED_NULL = -2, 340 KSP_DIVERGED_ITS = -3, 341 KSP_DIVERGED_DTOL = -4, 342 KSP_DIVERGED_BREAKDOWN = -5, 343 KSP_DIVERGED_BREAKDOWN_BICG = -6, 344 KSP_DIVERGED_NONSYMMETRIC = -7, 345 KSP_DIVERGED_INDEFINITE_PC = -8, 346 347 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 348 349 EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,int,PetscReal,KSPConvergedReason*,void*),void *); 350 EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 351 EXTERN PetscErrorCode KSPDefaultConverged(KSP,int,PetscReal,KSPConvergedReason*,void *); 352 EXTERN PetscErrorCode KSPSkipConverged(KSP,int,PetscReal,KSPConvergedReason*,void *); 353 EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 354 355 EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 356 357 /*E 358 KSPCGType - Determines what type of CG to use 359 360 Level: beginner 361 362 .seealso: KSPCGSetType() 363 E*/ 364 typedef enum {KSP_CG_SYMMETRIC=1,KSP_CG_HERMITIAN=2} KSPCGType; 365 366 EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 367 368 EXTERN PetscErrorCode PCPreSolve(PC,KSP); 369 EXTERN PetscErrorCode PCPostSolve(PC,KSP); 370 371 EXTERN PetscErrorCode KSPLGMonitorCreate(const char[],const char[],int,int,int,int,PetscDrawLG*); 372 EXTERN PetscErrorCode KSPLGMonitor(KSP,int,PetscReal,void*); 373 EXTERN PetscErrorCode KSPLGMonitorDestroy(PetscDrawLG); 374 EXTERN PetscErrorCode KSPLGTrueMonitorCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 375 EXTERN PetscErrorCode KSPLGTrueMonitor(KSP,int,PetscReal,void*); 376 EXTERN PetscErrorCode KSPLGTrueMonitorDestroy(PetscDrawLG); 377 378 PETSC_EXTERN_CXX_END 379 #endif 380