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