SVDTRLanczosSetOneSide#
Indicate if the variant of the Lanczos method to be used is one-sided or two-sided.
Synopsis#
#include "slepcsvd.h"
PetscErrorCode SVDTRLanczosSetOneSide(SVD svd,PetscBool oneside)
Logically Collective
Input Parameters#
svd - the singular value solver context
oneside - boolean flag indicating if the method is one-sided or not
Options Database Key#
-svd_trlanczos_oneside - enable the one-sided variant
Notes#
By default, a two-sided variant is selected, which is sometimes slightly more robust. However, the one-sided variant is faster because it avoids the orthogonalization associated to left singular vectors. See more details in [Hernandez et al., 2007, Hernández et al., 2008].
One-sided orthogonalization is also available for the GSVD, in which case two orthogonalizations out of three are avoided, see [Alvarruiz et al., 2024].
References#
F. Alvarruiz, C. Campos, and J. E. Roman. Thick-restarted joint Lanczos bidiagonalization for the GSVD. J. Comput. Appl. Math., 440:115506, 2024. doi:10.1016/j.cam.2023.115506.
V. Hernandez, J. E. Roman, and A. Tomas. Restarted Lanczos bidiagonalization for the SVD in SLEPc. Technical Report STR-8, Universitat Politècnica de València, 2007. URL: https://slepc.upv.es/documentation.
V. Hernández, J. E. Román, and A. Tomás. A robust and efficient parallel SVD solver based on restarted Lanczos bidiagonalization. Electron. Trans. Numer. Anal., 31:68–85, 2008. URL: https://etna.ricam.oeaw.ac.at/volumes/2001-2010/vol31/abstract.php?pages=68-85.
See Also#
SVD: Singular Value Decomposition, SVDTRLANCZOS, SVDLanczosSetOneSide()
Level#
advanced
Location#
Implementations#
SVDTRLanczosSetOneSide_TRLanczos() in src/svd/impls/trlanczos/trlanczos.c
Index of all SVD routines Table of Contents for all manual pages Index of all manual pages