SVDCrossSetExplicitMatrix#
Indicate if the eigensolver operator \(A^*A\) must be computed explicitly.
Synopsis#
#include "slepcsvd.h"
PetscErrorCode SVDCrossSetExplicitMatrix(SVD svd,PetscBool explicitmat)
Logically Collective
Input Parameters#
svd - the singular value solver context
explicitmat -
PETSC_TRUEif \(A^*A\) must be built explicitly
Options Database Key#
-svd_cross_explicitmatrix - toggle the explicit construction of the matrix
Notes#
In GSVD there are two cross product matrices, \(A^*A\) and \(B^*B\). In HSVD the expression for the cross product matrix is different, \(A^*\Omega A\). See Mathematical Background for details.
By default the matrices are not built explicitly, but handled as shell matrices,
see MATSHELL.
See Also#
SVD: Singular Value Decomposition, Mathematical Background, SVDCROSS, SVDCrossGetExplicitMatrix(), MATSHELL
Level#
advanced
Location#
Implementations#
SVDCrossSetExplicitMatrix_Cross() in src/svd/impls/cross/cross.c
Index of all SVD routines Table of Contents for all manual pages Index of all manual pages