KSPDGMRES#
Implements the deflated GMRES as defined in [EBP96] and [WP13]
Options Database Keys#
GMRES Options (inherited)#
-ksp_gmres_restart
- the number of Krylov directions to orthogonalize against-ksp_gmres_haptol
- sets the tolerance for “happy ending” (exact convergence)-ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of vectors are allocated as needed)
-ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
-ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
-ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the stability of the classical Gram-Schmidt orthogonalization.
-ksp_gmres_krylov_monitor - plot the Krylov space generated
DGMRES Options Database Keys#
-ksp_dgmres_eigen
- number of smallest eigenvalues to extract at each restart-ksp_dgmres_max_eigen <max_neig> - maximum number of eigenvalues that can be extracted during the iterative process
-ksp_dgmres_force - use the deflation at each restart; switch off the adaptive strategy.
-ksp_dgmres_view_deflation_vecs
- View the deflation vectors, where viewerspec is a key that can be parsed byPetscOptionsCreateViewer()
. If neig > 1, viewerspec should end with “:append”. No vectors will be viewed if the adaptive strategy chooses not to deflate, so -ksp_dgmres_force should also be given. The deflation vectors span a subspace that may be a good approximation of the subspace of smallest eigenvectors of the preconditioned operator, so this option can aid in understanding the performance of a preconditioner.
Notes#
Left and right preconditioning are supported, but not symmetric preconditioning. Complex arithmetic is not supported
In this implementation, the adaptive strategy allows switching to deflated GMRES when the stagnation occurs.
Contributed by#
Desire NUENTSA WAKAM, INRIA
References#
Jocelyne Erhel, Kevin Burrage, and Bert Pohl. Restarted GMRES preconditioned by deflation. Journal of computational and applied mathematics, 69(2):303–318, 1996.
Désiré Nuentsa Wakam and François Pacull. Memory efficient hybrid algebraic solvers for linear systems arising from compressible flows. Computers & Fluids, 80:158–167, 2013.
See Also#
KSP: Linear System Solvers, KSPCreate()
, KSPSetType()
, KSPType
, KSP
, KSPFGMRES
, KSPLGMRES
,
KSPGMRESSetRestart()
, KSPGMRESSetHapTol()
, KSPGMRESSetPreAllocateVectors()
, KSPGMRESSetOrthogonalization()
, KSPGMRESGetOrthogonalization()
,
KSPGMRESClassicalGramSchmidtOrthogonalization()
, KSPGMRESModifiedGramSchmidtOrthogonalization()
,
KSPGMRESCGSRefinementType
, KSPGMRESSetCGSRefinementType()
, KSPGMRESGetCGSRefinementType()
, KSPGMRESMonitorKrylov()
, KSPSetPCSide()
Level#
beginner
Location#
src/ksp/ksp/impls/gmres/dgmres/dgmres.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages