petsc-3.6.1 2015-08-06
Report Typos and Errors

KSPDGMRES

Implements the deflated GMRES as defined in [1,2]. In this implementation, the adaptive strategy allows to switch to the deflated GMRES when the stagnation occurs.

Options Database Keys

GMRES Options (inherited)

-ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
-ksp_gmres_haptol <tol> - 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 <never,ifneeded,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 <neig> - 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 <viewerspec> - View the deflation vectors, where viewerspec is a key that can be parsed by PetscOptionsGetViewer(). 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 yet supported

References

[1]Restarted GMRES preconditioned by deflation,J. Computational and Applied Mathematics, 69(1996), 303-318. [2]On the performance of various adaptive preconditioned GMRES strategies, 5(1998), 101-121. [3] D. NUENTSA WAKAM and F. PACULL, Memory Efficient Hybrid Algebraic Solvers for Linear Systems Arising from Compressible Flows, Computers and Fluids, In Press, http://dx.doi.org/10.1016/j.compfluid.2012.03.023

Contributed by: Desire NUENTSA WAKAM,INRIA

.seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), 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