KSPCGNE#
Applies the preconditioned conjugate gradient method to the normal equations without explicitly forming \(A^T*A\)
Options Database Key#
-ksp_cg_type <Hermitian or symmetric - (for complex matrices only) indicates the matrix is Hermitian or symmetric
Notes#
Eigenvalue computation routines including KSPSetComputeEigenvalues()
and KSPComputeEigenvalues()
will return information about the
spectrum of \(A^T*A\), rather than \(A\).
KSPCGNE
is a general-purpose non-symmetric method. It works well when the singular values are much better behaved than
eigenvalues. A unitary matrix is a classic example where KSPCGNE
converges in one iteration, but KSPGMRES
and KSPCGS
need N
iterations, see [NRT90]. If you intend to solve least squares problems, use KSPLSQR
.
This is NOT a different algorithm than used with KSPCG
, it merely uses that algorithm with the
matrix defined by \(A^T*A\) and preconditioner defined by \(B^T*B\) where \(B\) is the preconditioner for \(A\).
This method requires that one be able to apply the transpose of the preconditioner and operator as well as the operator and preconditioner. If the transpose of the preconditioner is not available then the preconditioner is used in its place so one ends up preconditioning \(A^T*A\) with \(B*B\). Seems odd?
This only supports left preconditioning.
Developer Note#
This object is subclassed off of KSPCG
, see the source code in src/ksp/ksp/impls/cg for comments on the structure of the code
References#
- NRT90
Noël M. Nachtigal, Satish C. Reddy, and Lloyd N. Trefethen. How fast are nonsymmetric matrix iterations? Technical Report 90-2, Massachusetts Institute of Technology, March 1990.
See Also#
KSP: Linear System Solvers, KSPCreate()
, KSPSetType()
, KSPType
, KSP
, ‘KSPCG’, KSPLSQR', 'KSPCGLS
,
KSPCGSetType()
, KSPBICG
, KSPSetComputeEigenvalues()
, KSPComputeEigenvalues()
Level#
beginner
Location#
src/ksp/ksp/impls/cg/cgne/cgne.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages