HLIBpro 3.1
Loading...
Searching...
No Matches
BLAS

Classes

class  Matrix< T_value >
 Standard dense matrix in basic linear algebra, i.e. BLAS/LAPACK. More...
 
class  TransposeView< T_matrix >
 Provide transposed view of a matrix. More...
 
class  AdjoinView< T_matrix >
 Provide adjoint view, e.g. conjugate transposed of a given matrix. More...
 
class  MatrixView< T_matrix >
 Provide generic view to a matrix, e.g. transposed or adjoint. More...
 
class  MatrixBase< T_derived >
 defines basic interface for matrices More...
 
class  MemBlock< T_value >
 Defines a reference countable memory block. More...
 
class  Range
 indexset with modified ctors More...
 
class  Vector< T_value >
 Standard vector in basic linear algebra, i.e. BLAS/LAPACK. More...
 

Functions

matop_t conjugate (const matop_t op)
 return conjugate of given matrix operation
 
matop_t transposed (const matop_t op)
 return transposed of given matrix operation
 
matop_t adjoint (const matop_t op)
 return adjoint of given matrix operation
 
matop_t apply_op (const matop_t op0, const matop_t op1)
 return operator corresponding to op0( op1 )
 

Vector Algebra

template<typename T1 , typename T2 >
std::enable_if_t< is_vector< T2 >::value &&std::is_same< T1, typename T2::value_t >::value, voidfill (const T1 f, T2 &x)
 fill vector with constant
 
template<typename T1 >
std::enable_if_t< is_vector< T1 >::value, value_type_t< T1 > > sum (T1 &x)
 fill vector with repeated function evaluation, x_i = f()
 
template<typename T1 >
std::enable_if_t< is_vector< T1 >::value, voidconj (T1 &x)
 conjugate entries in vector
 
template<typename T1 , typename T2 >
std::enable_if_t< is_vector< T2 >::value &&is_same_type< T1, typename T2::value_t >::value, voidscale (const T1 f, T2 &x)
 scale vector by constant
 
template<typename T1 , typename T2 >
std::enable_if_t< is_vector< T1 >::value &&is_vector< T2 >::value &&is_same_type< typename T1::value_t, typename T2::value_t >::value, voidcopy (const T1 &x, T2 &y)
 copy x into y
 
template<typename T1 , typename T2 >
std::enable_if_t< is_vector< T1 >::value &&is_vector< T2 >::value &&is_same_type< typename T1::value_t, typename T2::value_t >::value, voidswap (T1 &x, T2 &y)
 exchange x and y
 
template<typename T1 >
std::enable_if_t< is_vector< T1 >::value, idx_t > max_idx (const T1 &x)
 determine index with maximal absolute value in x
 
template<typename T1 >
std::enable_if_t< is_vector< T1 >::value, idx_t > min_idx (const T1 &x)
 determine index with minimax absolute value in x
 
template<typename T1 , typename T2 , typename T3 >
std::enable_if_t< is_vector< T2 >::value &&is_vector< T3 >::value &&is_same_type< T1, typename T2::value_t >::value &&is_same_type< T1, typename T3::value_t >::value, voidadd (const T1 alpha, const T2 &x, T3 &y)
 compute y ≔ y + α·x
 
template<typename T1 , typename T2 >
std::enable_if_t< is_vector< T1 >::value &&is_vector< T2 >::value &&is_same_type< typename T1::value_t, typename T2::value_t >::value, typename T1::value_t > dot (const T1 &x, const T2 &y)
 compute <x,y> = x^H · y
 
template<typename T1 , typename T2 >
std::enable_if_t< is_vector< T1 >::value &&is_vector< T2 >::value &&is_same_type< typename T1::value_t, typename T2::value_t >::value, typename T1::value_t > dotu (const T1 &x, const T2 &y)
 compute <x,y> without conjugating x, e.g. x^T · y
 
template<typename T1 >
std::enable_if< is_vector< T1 >::value, real_type_t< typenameT1::value_t > >::type norm2 (const T1 &x)
 compute ∥x∥₂
 
template<typename T1 , typename T2 >
std::enable_if_t< is_vector< T1 >::value &&is_vector< T2 >::value &&is_same_type< typename T1::value_t, typename T2::value_t >::value, typename T1::value_t > stable_dotu (const T1 &x, const T2 &y)
 compute dot product x · y numerically stable
 
template<typename T1 >
std::enable_if_t< is_vector< T1 >::value, typename T1::value_t > stable_sum (const T1 &x)
 compute sum of elements in x numerically stable
 

Basic Matrix Algebra

template<typename T1 >
void fill_rand (Matrix< T1 > &M)
 fill M with random entries
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidtranspose (T1 &A)
 transpose matrix A: A → A^T ASSUMPTION: A is square matrix
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidconj_transpose (T1 &A)
 conjugate transpose matrix A: A → A^H ASSUMPTION: A is square matrix
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidmax_idx (const T1 &M, idx_t &row, idx_t &col)
 determine index (i,j) with maximal absolute value in M and return in row and col
 
template<typename T1 , typename T2 , typename T3 , typename T4 >
std::enable_if_t< is_vector< T2 >::value &&is_vector< T3 >::value &&is_matrix< T4 >::value &&is_same_type< T1, typename T2::value_t >::value &&is_same_type< T1, typename T3::value_t >::value &&is_same_type< T1, typename T4::value_t >::value, voidadd_r1 (const T1 alpha, const T2 &x, const T3 &y, T4 &A)
 compute A ≔ A + α·x·y^H
 
template<typename T1 , typename T2 , typename T3 , typename T4 >
std::enable_if_t< is_vector< T2 >::value &&is_vector< T3 >::value &&is_matrix< T4 >::value &&is_same_type< T1, typename T2::value_t >::value &&is_same_type< T1, typename T3::value_t >::value &&is_same_type< T1, typename T4::value_t >::value, voidadd_r1u (const T1 alpha, const T2 &x, const T3 &y, T4 &A)
 compute A ≔ A + α·x·y^T
 
template<typename T1 , typename T2 , typename T3 , typename T4 >
std::enable_if_t< is_matrix< T2 >::value &&is_vector< T3 >::value &&is_vector< T4 >::value &&is_same_type< T1, typename T2::value_t >::value &&is_same_type< T1, typename T3::value_t >::value &&is_same_type< T1, typename T4::value_t >::value, voidmulvec (const T1 alpha, const T2 &A, const T3 &x, const T1 beta, T4 &y)
 compute y ≔ β·y + α·A·x
 
template<typename T1 , typename T2 , typename T3 >
std::enable_if_t< is_matrix< T2 >::value &&is_vector< T3 >::value &&is_same_type< T1, typename T2::value_t >::value &&is_same_type< T1, typename T3::value_t >::value, Vector< typename T2::value_t > > mulvec (const T1 alpha, const T2 &A, const T3 &x)
 compute y ≔ α·A·x
 
template<typename T1 , typename T2 >
std::enable_if_t< is_matrix< T1 >::value &&is_vector< T2 >::value &&is_same_type< typename T1::value_t, typename T2::value_t >::value, voidmulvec_tri (const tri_type_t shape, const diag_type_t diag, const T1 &A, T2 &x)
 compute x ≔ M·x, where M is upper or lower triangular with unit or non-unit diagonal
 
template<typename T1 , typename T2 >
std::enable_if_t< is_matrix< T1 >::value &&is_vector< T2 >::value &&is_same_type< typename T1::value_t, typename T2::value_t >::value, voidsolve (T1 &A, T2 &b)
 solve A·x = b with known A and b; x overwrites b (A is overwritten upon exit!)
 
template<typename T1 , typename T2 , typename T3 , typename T4 >
std::enable_if_t< is_matrix< T2 >::value &&is_matrix< T3 >::value &&is_matrix< T4 >::value &&is_same_type< T1, typename T2::value_t >::value &&is_same_type< T1, typename T3::value_t >::value &&is_same_type< T1, typename T4::value_t >::value, voidprod (const T1 alpha, const T2 &A, const T3 &B, const T1 beta, T4 &C)
 compute C ≔ β·C + α·A·B
 
template<typename T1 , typename T2 , typename T3 >
std::enable_if_t< is_matrix< T2 >::value &&is_matrix< T3 >::value &&is_same_type< T1, typename T2::value_t >::value &&is_same_type< T1, typename T3::value_t >::value, Matrix< typename T2::value_t > > prod (const T1 alpha, const T2 &A, const T3 &B)
 compute C ≔ α·A·B
 
template<typename T1 , typename T2 >
std::enable_if_t< is_matrix< T1 >::value &&is_matrix< T2 >::value &&is_same_type< typename T1::value_t, typename T2::value_t >::value, voidhadamard_prod (const T1 &A, T2 &B)
 compute B ≔ A⊙B, e.g. Hadamard product
 
template<typename T1 , typename T2 , typename T3 >
std::enable_if_t< is_matrix< T2 >::value &&is_matrix< T3 >::value &&is_same_type< T1, typename T2::value_t >::value &&is_same_type< T1, typename T3::value_t >::value, voidprod_tri (const eval_side_t side, const tri_type_t uplo, const diag_type_t diag, const T1 alpha, const T2 &A, T3 &B)
 compute B ≔ α·A·B or B ≔ α·B·A with triangular matrix A
 
template<typename T1 , typename T2 >
std::enable_if_t< is_matrix< T1 >::value &&is_vector< T2 >::value &&is_same_type< real_type_t< typename T1::value_t >, real_type_t< typename T2::value_t > >::value, voidprod_diag (T1 &M, const T2 &D, const idx_t k)
 multiply k columns of M with diagonal matrix D, e.g. compute M ≔ M·D
 
template<typename T1 , typename T2 >
std::enable_if_t< is_vector< T1 >::value &&is_matrix< T2 >::value &&is_same_type< real_type_t< typename T1::value_t >, real_type_t< typename T2::value_t > >::value, voidprod_diag (const T1 &D, T2 &M, const idx_t k)
 multiply k rows of M with diagonal matrix D, e.g. compute M ≔ D·M
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, real_type_t< typename T1::value_t > > normF (const T1 &M)
 return Frobenius norm of M
 
template<typename T1 , typename T2 >
std::enable_if_t< is_matrix< T1 >::value &&is_matrix< T2 >::value &&is_same_type< typename T1::value_t, typename T2::value_t >::value, real_type_t< typename T1::value_t > > diff_normF (const T1 &A, const T2 &B)
 compute Frobenius norm of A-B
 
template<typename T1 , typename T2 >
std::enable_if_t< is_matrix< T1 >::value &&is_matrix< T2 >::value &&is_same_type< typename T1::value_t, typename T2::value_t >::value, real_type_t< typename T1::value_t > > lr_normF (const T1 &A, const T2 &B)
 compute Frobenius norm of M=A·B^H
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, real_type_t< typename T1::value_t > > cond (const T1 &M)
 return condition of M
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidmake_symmetric (T1 &A)
 make given matrix symmetric, e.g. copy lower left part to upper right part
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidmake_hermitian (T1 &A)
 make given matrix hermitian, e.g. copy conjugated lower left part to upper right part and make diagonal real
 

Advanced Matrix Algebra

template<typename T1 , typename T2 >
std::enable_if_t< is_matrix< T1 >::value &&is_vector< T2 >::value &&is_same_type< typename T1::value_t, typename T2::value_t >::value, voidsolve_tri (const tri_type_t uplo, const diag_type_t diag, const T1 &A, T2 &b)
 solve A·x = b with known A and b; x overwrites b
 
template<typename T1 , typename T2 , typename T3 >
std::enable_if_t< is_matrix< T2 >::value &&is_matrix< T3 >::value &&is_same_type< T1, typename T2::value_t >::value &&is_same_type< T1, typename T3::value_t >::value, voidsolve_tri (const eval_side_t side, const tri_type_t uplo, const diag_type_t diag, const T1 alpha, const T2 &A, T3 &B)
 solve A·X = α·B or X·A· = α·B with known A and B; X overwrites B
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidinvert (T1 &A)
 invert matrix A; A will be overwritten with A^-1 upon exit
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidinvert (T1 &A, const tri_type_t tri_type, const diag_type_t diag_type)
 invert lower or upper triangular matrix A with unit or non-unit diagonal; A will be overwritten with A^-1 upon exit
 
template<typename T >
void pseudo_invert (Matrix< T > &A, const TTruncAcc &acc)
 compute pseudo inverse of matrix A with precision acc
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidlu (T1 &A)
 compute LU factorisation of the n×m matrix A with n×min(n,m) unit diagonal lower triangular matrix L and min(n,m)xm upper triangular matrix U; A will be overwritten with L and U upon exit
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidllt (T1 &A)
 compute L·L^T factorisation of given symmetric n×n matrix A
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidllh (T1 &A)
 compute L·L^H factorisation of given hermitian n×n matrix A
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidldlt (T1 &A)
 compute L·D·L^T factorisation of given symmetric n×n matrix A
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidldlh (T1 &A)
 compute L·D·L^H factorisation of given hermitian n×n matrix A
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidqr (T1 &A, Matrix< typename T1::value_t > &R)
 Compute QR factorisation of the matrix A.
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidtsqr (T1 &A, Matrix< typename T1::value_t > &R, const size_t ntile=0)
 Compute QR factorisation of the tall-and-skinny matrix A.
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidqrp (T1 &A, Matrix< typename T1::value_t > &R, std::vector< blas_int_t > &P)
 Compute QR factorisation with column pivoting of the matrix A.
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voideigen (T1 &M, Vector< typename T1::value_t > &eig_val, Matrix< typename T1::value_t > &eig_vec)
 compute eigenvalues and eigenvectors of matrix M
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voideigen_herm (T1 &M, Vector< real_type_t< typename T1::value_t > > &eig_val, Matrix< typename T1::value_t > &eig_vec)
 compute eigenvalues and eigenvectors of the hermitian matrix M
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voideigen (T1 &M, const Range &eig_range, Vector< typename T1::value_t > &eig_val, Matrix< typename T1::value_t > &eig_vec)
 compute selected (by eig_range) eigenvalues and eigenvectors of the symmetric matrix M
 
template<typename T1 , typename T2 >
std::enable_if_t< is_vector< T1 >::value &&is_vector< T2 >::value &&is_same_type< typename T1::value_t, typename T2::value_t >::value, voideigen (T1 &diag, T2 &subdiag, Vector< typename T1::value_t > &eig_val, Matrix< typename T1::value_t > &eig_vec)
 compute eigenvalues and eigenvectors of the symmetric, tridiagonal matrix defines by diagonal coefficients in diag and off-diagonal coefficients subdiag
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidsvd (T1 &A, Vector< real_type_t< typename T1::value_t > > &S, Matrix< typename T1::value_t > &V)
 compute SVD decomposition \( A = U·S·V^H \) of the nxm matrix A with n×min(n,m) matrix U, min(n,m)×min(n,m) matrix S (diagonal) and m×min(n,m) matrix V; A will be overwritten with U upon exit
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidsvd (T1 &A, Vector< real_type_t< typename T1::value_t > > &S, const bool left=true)
 compute SVD decomposition \( A = U·S·V^H \) of the nxm matrix A but return only the left/right singular vectors and the singular values S ∈ ℝ^min(n,m); upon exit, A will be contain the corresponding sing. vectors
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidsv (T1 &A, Vector< real_type_t< typename T1::value_t > > &S)
 compute SVD decomposition \( A = U·S·V^H \) of the nxm matrix A but return only the singular values S ∈ ℝ^min(n,m); A will be overwritten with U upon exit
 
template<typename T1 , typename T2 >
std::enable_if_t< is_matrix< T1 >::value &&is_matrix< T2 >::value &&is_same_type< typename T1::value_t, typename T2::value_t >::value, voidsv (T1 &A, T2 &B, Vector< real_type_t< typename T1::value_t > > &S)
 compute SVD decomposition \( M = A·B^H = U·S·V^H \) of the nxm low-rank matrix M but return only the singular values S ∈ ℝ^min(n,m); A and B will be overwritten upon exit
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, size_tapprox (T1 &M, const TTruncAcc &acc, Matrix< typename T1::value_t > &A, Matrix< typename T1::value_t > &B)
 approximate given dense matrix M by low rank matrix according to accuracy acc. The low rank matrix will be stored in A and B
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, size_tapprox_svd (T1 &M, const TTruncAcc &acc, Matrix< typename T1::value_t > &A, Matrix< typename T1::value_t > &B)
 approximate given dense matrix M by low rank matrix according to accuracy acc using SVD. The low rank matrix will be stored in A and B
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, size_tapprox_rrqr (T1 &M, const TTruncAcc &acc, Matrix< typename T1::value_t > &A, Matrix< typename T1::value_t > &B)
 approximate given dense matrix M by low rank matrix according to accuracy acc using RRQR. The low rank matrix will be stored in A and B.
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, size_tapprox_randsvd (T1 &M, const TTruncAcc &acc, Matrix< typename T1::value_t > &A, Matrix< typename T1::value_t > &B, const uint power_steps=CFG::BLAS::power_steps, const uint oversampling=CFG::BLAS::oversampling)
 approximate given dense matrix M by low rank matrix according to accuracy acc using randomized SVD. The low rank matrix will be stored in A and B.
 
template<typename T >
size_t truncate (Matrix< T > &A, Matrix< T > &B, const TTruncAcc &acc)
 truncate A · B^H based on accuracy acc.
 
template<typename T >
size_t truncate_svd (Matrix< T > &A, Matrix< T > &B, const TTruncAcc &acc)
 truncate A · B^H based on accuracy acc using SVD.
 
template<typename T >
size_t truncate_rrqr (Matrix< T > &A, Matrix< T > &B, const TTruncAcc &acc)
 truncate A · B^H based on accuracy acc using RRQR.
 
template<typename T >
size_t truncate_rand (Matrix< T > &A, Matrix< T > &B, const TTruncAcc &acc)
 truncate A · B^H based on accuracy acc using randomized SVD.
 
template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, voidfactorise_ortho (T1 &A, Matrix< typename T1::value_t > &R)
 construct factorisation A = Q·R of A, with orthonormal Q
 
template<typename T >
void factorise_ortho (Matrix< T > &A, Matrix< T > &R, const TTruncAcc &acc)
 construct approximate factorisation A = Q·R of A, with orthonormal Q
 
void print_statistics ()
 print statistics for Algebra functions
 
void reset_statistics ()
 reset statistics for Algebra functions
 

Matrix Modifiers

Classes for matrix modifiers, e.g. transposed and adjoint view.

template<typename T >
TransposeView< T > transposed (const T &M)
 return transposed view object for matrix
 
template<typename T >
AdjoinView< T > adjoint (const T &M)
 return adjoint view object for matrix
 
template<typename T >
MatrixView< T > mat_view (const matop_t op, const T &M)
 convert matop_t into view object
 

Detailed Description

This modules provides most low level algebra functions, e.g. vector dot products, matrix multiplication, factorisation and singular value decomposition. See also BLAS/LAPACK Interface for an introduction.

To include all BLAS algebra functions and classes add

#include <hlib-blas.hh>

to your source files.

Function Documentation

◆ adjoint() [1/2]

matop_t adjoint ( const matop_t  op)
inline
Parameters
opmatrix op. to be adjoint

◆ adjoint() [2/2]

template<typename T >
AdjoinView< T > adjoint ( const T &  M)
Parameters
Mmatrix to conjugate transpose

◆ cond()

template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, real_type_t< typename T1::value_t > > cond ( const T1 M)
  • if M ≡ 0 or size(M) = 0, 0 is returned

◆ conjugate()

matop_t conjugate ( const matop_t  op)
inline
Parameters
opmatrix op. to be conjugated

◆ eigen()

template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, void > eigen ( T1 M,
const Range eig_range,
Vector< typename T1::value_t > &  eig_val,
Matrix< typename T1::value_t > &  eig_vec 
)
  • the lower half of M is accessed

◆ factorise_ortho() [1/2]

template<typename T >
void factorise_ortho ( Matrix< T > &  A,
Matrix< T > &  R,
const TTruncAcc acc 
)
  • approximation quality is definded by acc
  • A is overwritten with Q upon exit

◆ factorise_ortho() [2/2]

template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, void > factorise_ortho ( T1 A,
Matrix< typename T1::value_t > &  R 
)
  • A is overwritten with Q upon exit

◆ mat_view()

template<typename T >
MatrixView< T > mat_view ( const matop_t  op,
const T &  M 
)
Parameters
opmatop_t value
Mmatrix to create view for

◆ norm2()

template<typename T1 >
std::enable_if< is_vector< T1 >::value, real_type_t< typenameT1::value_t > >::type norm2 ( const T1 x)

return spectral norm of M

◆ pseudo_invert()

template<typename T >
void pseudo_invert ( Matrix< T > &  A,
const TTruncAcc acc 
)

\detail Compute pseudo inverse B of matrix A up to precision acc, e.g. \(\|A-B\|\le \epsilon\) with \(\epsilon\) defined by acc.

◆ qr()

template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, void > qr ( T1 A,
Matrix< typename T1::value_t > &  R 
)
   Compute QR factorisation of the n×m matrix \a A with
   n×m matrix Q and mxm matrix R (n >= m); \a A will be
   overwritten with Q upon exit

◆ qrp()

template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, void > qrp ( T1 A,
Matrix< typename T1::value_t > &  R,
std::vector< blas_int_t > &  P 
)
   Compute QR factorisation with column pivoting of the n×m matrix \a A,
   i.e., \f$ AP = QR \f$, with n×m matrix Q and mxm matrix R (n >= m).
   \a P_i = j means, column i of AP was column j of A; \a A will be
   overwritten with Q upon exit

◆ stable_dotu()

template<typename T1 , typename T2 >
std::enable_if_t< is_vector< T1 >::value && is_vector< T2 >::value && is_same_type< typename T1::value_t, typename T2::value_t >::value, typename T1::value_t > stable_dotu ( const T1 x,
const T2 y 
)
Parameters
xfirst argument of dot product
ysecond argument of dot product

◆ stable_sum()

template<typename T1 >
std::enable_if_t< is_vector< T1 >::value, typename T1::value_t > stable_sum ( const T1 x)
Parameters
xvector holding coefficients to sum up

◆ sum()

template<typename T1 >
std::enable_if_t< is_vector< T1 >::value, value_type_t< T1 > > sum ( T1 x)

sum elements in vector

◆ transposed() [1/2]

matop_t transposed ( const matop_t  op)
inline
Parameters
opmatrix op. to be transposed

◆ transposed() [2/2]

template<typename T >
TransposeView< T > transposed ( const T &  M)
Parameters
Mmatrix to transpose

◆ truncate()

template<typename T >
size_t truncate ( Matrix< T > &  A,
Matrix< T > &  B,
const TTruncAcc acc 
)
   Truncate given \a A · \a B^H low rank matrix matrix 
   (\a A being n×k and \a B being m×k) with respect to
   given accuracy \a acc; store truncated matrix in
   A(:,0:k-1) and B(:,0:k-1) where k is the returned
   new rank after truncation

◆ tsqr()

template<typename T1 >
std::enable_if_t< is_matrix< T1 >::value, void > tsqr ( T1 A,
Matrix< typename T1::value_t > &  R,
const size_t  ntile = 0 
)
   Compute QR factorisation of the n×m matrix \a A, m ≪ n, with
   n×m matrix Q and mxm matrix R (n >= m); \a A will be
   overwritten with Q upon exit
   - ntile defines tile size (0: use internal default tile size)