![]() |
Blaze 3.9
|
Modules | |
Expressions | |
CompressedMatrix | |
IdentityMatrix | |
ZeroMatrix | |
Classes | |
class | blaze::SparseMatrix< MT, SO > |
Base class for sparse matrices. More... | |
class | blaze::MatrixAccessProxy< MT > |
Access proxy for sparse, ![]() | |
Functions | |
template<typename MT1 , bool SO1, typename MT2 , bool SO2> | |
decltype(auto) | blaze::kron (const DenseMatrix< MT1, SO1 > &lhs, const SparseMatrix< MT2, SO2 > &rhs) |
Computes the Kronecker product of a dense matrix and a sparse matrix ( ![]() | |
template<typename VT1 , typename VT2 > | |
decltype(auto) | blaze::operator* (const DenseVector< VT1, false > &lhs, const SparseVector< VT2, true > &rhs) |
Multiplication operator for the dense vector-sparse vector outer product ( ![]() | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::decldiag (const SparseMatrix< MT, SO > &sm) |
Declares the given sparse matrix expression sm as diagonal. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::declherm (const SparseMatrix< MT, SO > &sm) |
Declares the given sparse matrix expression sm as Hermitian. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::decllow (const SparseMatrix< MT, SO > &sm) |
Declares the given sparse matrix expression sm as lower. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::declstrlow (const SparseMatrix< MT, SO > &sm) |
Declares the given sparse matrix expression sm as strictly lower. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::declstrupp (const SparseMatrix< MT, SO > &sm) |
Declares the given sparse matrix expression sm as strictly upper. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::declsym (const SparseMatrix< MT, SO > &sm) |
Declares the given sparse matrix expression sm as symmetric. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::declunilow (const SparseMatrix< MT, SO > &sm) |
Declares the given sparse matrix expression sm as unilower. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::decluniupp (const SparseMatrix< MT, SO > &sm) |
Declares the given sparse matrix expression sm as uniupper. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::declupp (const SparseMatrix< MT, SO > &sm) |
Declares the given sparse matrix expression sm as upper. More... | |
template<typename MT1 , bool SO1, typename MT2 , bool SO2> | |
decltype(auto) | blaze::kron (const SparseMatrix< MT1, SO1 > &lhs, const DenseMatrix< MT2, SO2 > &rhs) |
Computes the Kronecker product of a sparse matrix and a dense matrix ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator% (const SparseMatrix< MT1, false > &lhs, const DenseMatrix< MT2, false > &rhs) |
Operator for the Schur product of a row-major sparse matrix and a row-major dense matrix ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator% (const SparseMatrix< MT1, false > &lhs, const DenseMatrix< MT2, true > &rhs) |
Operator for the Schur product of a row-major sparse matrix and a column-major dense matrix ( ![]() | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::eval (const SparseMatrix< MT, SO > &sm) |
Forces the evaluation of the given sparse matrix expression sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::fix (SparseMatrix< MT, SO > &sm) noexcept |
Fixing the size of the given sparse matrix. More... | |
template<typename MT , bool SO, typename OP > | |
decltype(auto) | blaze::map (const SparseMatrix< MT, SO > &sm, OP op) |
Evaluates the given custom operation on each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO, typename OP > | |
decltype(auto) | blaze::forEach (const SparseMatrix< MT, SO > &sm, OP op) |
Evaluates the given custom operation on each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::abs (const SparseMatrix< MT, SO > &sm) |
Applies the abs() function to each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::sign (const SparseMatrix< MT, SO > &sm) |
Applies the sign() function to each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::floor (const SparseMatrix< MT, SO > &sm) |
Applies the floor() function to each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::ceil (const SparseMatrix< MT, SO > &sm) |
Applies the ceil() function to each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::trunc (const SparseMatrix< MT, SO > &sm) |
Applies the trunc() function to each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::round (const SparseMatrix< MT, SO > &sm) |
Applies the round() function to each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::conj (const SparseMatrix< MT, SO > &sm) |
Returns a matrix containing the complex conjugate of each single element of sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::ctrans (const SparseMatrix< MT, SO > &sm) |
Returns the conjugate transpose matrix of sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::real (const SparseMatrix< MT, SO > &sm) |
Returns a matrix containing the real parts of each single element of sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::imag (const SparseMatrix< MT, SO > &sm) |
Returns a matrix containing the imaginary parts of each single element of sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::arg (const SparseMatrix< MT, SO > &sm) |
Returns a matrix containing the phase angle of each single element of sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::sqrt (const SparseMatrix< MT, SO > &sm) |
Computes the square root of each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::invsqrt (const SparseMatrix< MT, SO > &sm) |
Computes the inverse square root of each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::cbrt (const SparseMatrix< MT, SO > &sm) |
Computes the cubic root of each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::invcbrt (const SparseMatrix< MT, SO > &sm) |
Computes the inverse cubic root of each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO, typename DT > | |
decltype(auto) | blaze::clamp (const SparseMatrix< MT, SO > &sm, const DT &min, const DT &max) |
Restricts each single element of the sparse matrix sm to the range ![]() | |
template<typename MT , bool SO, typename ST , EnableIf_t< IsScalar_v< ST > > * = nullptr> | |
decltype(auto) | blaze::pow (const SparseMatrix< MT, SO > &sm, ST exp) |
Computes the exponential value for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::exp (const SparseMatrix< MT, SO > &sm) |
Computes ![]() | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::exp2 (const SparseMatrix< MT, SO > &sm) |
Computes ![]() | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::exp10 (const SparseMatrix< MT, SO > &sm) |
Computes ![]() | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::log (const SparseMatrix< MT, SO > &sm) |
Computes the natural logarithm for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::log2 (const SparseMatrix< MT, SO > &sm) |
Computes the binary logarithm for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::log10 (const SparseMatrix< MT, SO > &sm) |
Computes the common logarithm for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::log1p (const SparseMatrix< MT, SO > &sm) |
Computes the natural logarithm of x+1 for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::lgamma (const SparseMatrix< MT, SO > &sm) |
Computes the natural logarithm of the absolute value of the gamma function for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::sin (const SparseMatrix< MT, SO > &sm) |
Computes the sine for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::asin (const SparseMatrix< MT, SO > &sm) |
Computes the inverse sine for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::sinh (const SparseMatrix< MT, SO > &sm) |
Computes the hyperbolic sine for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::asinh (const SparseMatrix< MT, SO > &sm) |
Computes the inverse hyperbolic sine for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::cos (const SparseMatrix< MT, SO > &sm) |
Computes the cosine for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::acos (const SparseMatrix< MT, SO > &sm) |
Computes the inverse cosine for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::cosh (const SparseMatrix< MT, SO > &sm) |
Computes the hyperbolic cosine for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::acosh (const SparseMatrix< MT, SO > &sm) |
Computes the inverse hyperbolic cosine for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::tan (const SparseMatrix< MT, SO > &sm) |
Computes the tangent for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::atan (const SparseMatrix< MT, SO > &sm) |
Computes the inverse tangent for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::tanh (const SparseMatrix< MT, SO > &sm) |
Computes the hyperbolic tangent for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::atanh (const SparseMatrix< MT, SO > &sm) |
Computes the inverse hyperbolic tangent for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::erf (const SparseMatrix< MT, SO > &sm) |
Computes the error function for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::erfc (const SparseMatrix< MT, SO > &sm) |
Computes the complementary error function for each non-zero element of the sparse matrix sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::mean (const SparseMatrix< MT, SO > &sm) |
Computes the (arithmetic) mean for the given sparse matrix. More... | |
template<ReductionFlag , typename MT , bool SO> | |
decltype(auto) | blaze::mean (const SparseMatrix< MT, SO > &sm) |
Computes the row-/columnwise mean function for the given sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::noalias (const SparseMatrix< MT, SO > &sm) |
Forces the non-aliased evaluation of the given sparse matrix expression sm. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::norm (const SparseMatrix< MT, SO > &sm) |
Computes the L2 norm for the given sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::sqrNorm (const SparseMatrix< MT, SO > &sm) |
Computes the squared L2 norm for the given sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::l1Norm (const SparseMatrix< MT, SO > &sm) |
Computes the L1 norm for the given sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::l2Norm (const SparseMatrix< MT, SO > &sm) |
Computes the L2 norm for the given sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::l3Norm (const SparseMatrix< MT, SO > &sm) |
Computes the L3 norm for the given sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::l4Norm (const SparseMatrix< MT, SO > &sm) |
Computes the L4 norm for the given sparse matrix. More... | |
template<typename MT , bool SO, typename ST > | |
decltype(auto) | blaze::lpNorm (const SparseMatrix< MT, SO > &sm, ST p) |
Computes the Lp norm for the given sparse matrix. More... | |
template<size_t P, typename MT , bool SO> | |
decltype(auto) | blaze::lpNorm (const SparseMatrix< MT, SO > &sm) |
Computes the Lp norm for the given sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::linfNorm (const SparseMatrix< MT, SO > &sm) |
Computes the infinity norm for the given sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::maxNorm (const SparseMatrix< MT, SO > &sm) |
Computes the maximum norm for the given sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::minNorm (const SparseMatrix< MT, SO > &sm) |
Computes the minimum norm for the given sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::nosimd (const SparseMatrix< MT, SO > &sm) |
Disables the SIMD evaluation of the given sparse matrix expression sm. More... | |
template<typename MT , bool SO, typename OP > | |
decltype(auto) | blaze::reduce (const SparseMatrix< MT, SO > &sm, OP op) |
Performs a custom reduction operation on the given sparse matrix. More... | |
template<ReductionFlag , typename MT , bool SO, typename OP > | |
decltype(auto) | blaze::reduce (const SparseMatrix< MT, SO > &sm, OP op) |
Performs a custom reduction operation on the given sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::sum (const SparseMatrix< MT, SO > &sm) |
Reduces the given sparse matrix by means of addition. More... | |
template<ReductionFlag RF, typename MT , bool SO> | |
decltype(auto) | blaze::sum (const SparseMatrix< MT, SO > &sm) |
Reduces the given sparse matrix by means of addition. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::prod (const SparseMatrix< MT, SO > &sm) |
Reduces the given sparse matrix by means of multiplication. More... | |
template<ReductionFlag RF, typename MT , bool SO> | |
decltype(auto) | blaze::prod (const SparseMatrix< MT, SO > &sm) |
Reduces the given sparse matrix by means of multiplication. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::min (const SparseMatrix< MT, SO > &sm) |
Returns the smallest element of the sparse matrix. More... | |
template<ReductionFlag RF, typename MT , bool SO> | |
decltype(auto) | blaze::min (const SparseMatrix< MT, SO > &sm) |
Returns the smallest element of each row/columns of the sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::max (const SparseMatrix< MT, SO > &sm) |
Returns the largest element of the sparse matrix. More... | |
template<ReductionFlag RF, typename MT , bool SO> | |
decltype(auto) | blaze::max (const SparseMatrix< MT, SO > &sm) |
Returns the largest element of each row/columns of the sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::repeat (const SparseMatrix< MT, SO > &sm, size_t m, size_t n) |
Repeats the given sparse matrix. More... | |
template<size_t R0, size_t R1, typename MT , bool SO> | |
decltype(auto) | blaze::repeat (const SparseMatrix< MT, SO > &sm) |
Repeats the given sparse matrix. More... | |
template<typename MT , bool SO, typename ST , EnableIf_t< IsScalar_v< ST > > * = nullptr> | |
decltype(auto) | blaze::operator/ (const SparseMatrix< MT, SO > &mat, ST scalar) |
Division operator for the division of a sparse matrix by a scalar value ( ![]() | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::operator- (const SparseMatrix< MT, SO > &sm) |
Unary minus operator for the negation of a sparse matrix ( ![]() | |
template<typename MT , bool SO, typename ST , EnableIf_t< IsScalar_v< ST > > * = nullptr> | |
decltype(auto) | blaze::operator* (const SparseMatrix< MT, SO > &mat, ST scalar) |
Multiplication operator for the multiplication of a sparse matrix and a scalar value ( ![]() | |
template<typename ST , typename MT , bool SO, EnableIf_t< IsScalar_v< ST > > * = nullptr> | |
decltype(auto) | blaze::operator* (ST scalar, const SparseMatrix< MT, SO > &mat) |
Multiplication operator for the multiplication of a scalar value and a sparse matrix ( ![]() | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::serial (const SparseMatrix< MT, SO > &sm) |
Forces the serial evaluation of the given sparse matrix expression sm. More... | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator+ (const SparseMatrix< MT1, false > &lhs, const SparseMatrix< MT2, false > &rhs) |
Addition operator for the addition of two row-major sparse matrices ( ![]() | |
template<typename MT1 , bool SO1, typename MT2 , bool SO2> | |
bool | blaze::operator== (const SparseMatrix< MT1, SO1 > &lhs, const SparseMatrix< MT2, SO2 > &rhs) |
Equality operator for the comparison of two sparse matrices. More... | |
template<typename MT1 , bool SO1, typename MT2 , bool SO2> | |
bool | blaze::operator!= (const SparseMatrix< MT1, SO1 > &lhs, const SparseMatrix< MT2, SO2 > &rhs) |
Inequality operator for the comparison of two sparse matrices. More... | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::kron (const SparseMatrix< MT1, false > &lhs, const SparseMatrix< MT2, false > &rhs) |
Computes the Kronecker product of two row-major sparse matrices ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator* (const SparseMatrix< MT1, false > &lhs, const SparseMatrix< MT2, false > &rhs) |
Multiplication operator for the multiplication of two row-major sparse matrices ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator% (const SparseMatrix< MT1, false > &lhs, const SparseMatrix< MT2, false > &rhs) |
Operator for the Schur product of two row-major sparse matrices ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator- (const SparseMatrix< MT1, false > &lhs, const SparseMatrix< MT2, false > &rhs) |
Subtraction operator for the subtraction of two row-major sparse matrices ( ![]() | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::stddev (const SparseMatrix< MT, SO > &sm) |
Computes the standard deviation for the given sparse matrix. More... | |
template<ReductionFlag , typename MT , bool SO> | |
decltype(auto) | blaze::stddev (const SparseMatrix< MT, SO > &sm) |
Computes the row-/columnwise standard deviation function for the given sparse matrix. More... | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::trans (const SparseMatrix< MT, SO > &sm) |
Calculation of the transpose of the given sparse matrix. More... | |
template<bool B, typename MT , bool SO> | |
decltype(auto) | blaze::transIf (const SparseMatrix< MT, SO > &sm) |
Conditional calculation of the transpose of the given sparse matrix. More... | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator+ (const SparseMatrix< MT1, false > &lhs, const SparseMatrix< MT2, true > &rhs) |
Addition operator for the addition of a row-major and a column-major sparse matrix ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator+ (const SparseMatrix< MT1, true > &lhs, const SparseMatrix< MT2, false > &rhs) |
Addition operator for the addition of a column-major and a row-major sparse matrix ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::kron (const SparseMatrix< MT1, false > &lhs, const SparseMatrix< MT2, true > &rhs) |
Operator for the Kronecker product of a row-major and a column-major sparse matrix ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator* (const SparseMatrix< MT1, false > &lhs, const SparseMatrix< MT2, true > &rhs) |
Multiplication operator for the multiplication of a row-major sparse matrix and a column-major sparse matrix ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator% (const SparseMatrix< MT1, false > &lhs, const SparseMatrix< MT2, true > &rhs) |
Operator for the Schur product of a row-major and a column-major sparse matrix ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator- (const SparseMatrix< MT1, false > &lhs, const SparseMatrix< MT2, true > &rhs) |
Subtraction operator for the subtraction of a row-major and a column-major sparse matrix ( ![]() | |
template<typename MT , bool SO> | |
decltype(auto) | blaze::var (const SparseMatrix< MT, SO > &sm) |
Computes the variance for the given sparse matrix. More... | |
template<ReductionFlag , typename MT , bool SO> | |
decltype(auto) | blaze::var (const SparseMatrix< MT, SO > &sm) |
Computes the row-/column-wise variance function for the given sparse matrix. More... | |
template<typename VT1 , typename VT2 > | |
decltype(auto) | blaze::operator* (const SparseVector< VT1, false > &lhs, const DenseVector< VT2, true > &rhs) |
Multiplication operator for the sparse vector-dense vector outer product ( ![]() | |
template<typename VT1 , typename VT2 > | |
decltype(auto) | blaze::operator* (const SparseVector< VT1, false > &lhs, const SparseVector< VT2, true > &rhs) |
Multiplication operator for the sparse vector-sparse vector outer product ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator% (const SparseMatrix< MT1, true > &lhs, const DenseMatrix< MT2, false > &rhs) |
Operator for the Schur product of a column-major sparse matrix and a row-major dense matrix ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator% (const SparseMatrix< MT1, true > &lhs, const DenseMatrix< MT2, true > &rhs) |
Operator for the Schur product of a column-major sparse matrix and a column-major dense matrix ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::kron (const SparseMatrix< MT1, true > &lhs, const SparseMatrix< MT2, false > &rhs) |
Operator for the Kronecker product of a column-major and a row-major sparse matrix ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator* (const SparseMatrix< MT1, true > &lhs, const SparseMatrix< MT2, false > &rhs) |
Multiplication operator for the multiplication of a column-major sparse matrix and a row-major sparse matrix ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator% (const SparseMatrix< MT1, true > &lhs, const SparseMatrix< MT2, false > &rhs) |
Operator for the Schur product of a column-major and a row-major sparse matrix ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator- (const SparseMatrix< MT1, true > &lhs, const SparseMatrix< MT2, false > &rhs) |
Subtraction operator for the subtraction of a column-major and a row-major sparse matrix ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator+ (const SparseMatrix< MT1, true > &lhs, const SparseMatrix< MT2, true > &rhs) |
Addition operator for the addition of two column-major sparse matrices ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::kron (const SparseMatrix< MT1, true > &lhs, const SparseMatrix< MT2, true > &rhs) |
Computes the Kronecker product of two column-major sparse matrices ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator* (const SparseMatrix< MT1, true > &lhs, const SparseMatrix< MT2, true > &rhs) |
Multiplication operator for the multiplication of two column-major sparse matrices ( ![]() | |
template<typename MT1 , typename MT2 , DisableIf_t<(IsUniLower_v< MT1 > &&IsUniUpper_v< MT2 >)||(IsUniUpper_v< MT1 > &&IsUniLower_v< MT2 >)||(IsStrictlyLower_v< MT1 > &&IsUpper_v< MT2 >)||(IsStrictlyUpper_v< MT1 > &&IsLower_v< MT2 >)||(IsLower_v< MT1 > &&IsStrictlyUpper_v< MT2 >)||(IsUpper_v< MT1 > &&IsStrictlyLower_v< MT2 >)||(IsZero_v< MT1 >||IsZero_v< MT2 >) > * = nullptr> | |
const TSMatTSMatSchurExpr< MT1, MT2 > | blaze::tsmattsmatschur (const SparseMatrix< MT1, true > &lhs, const SparseMatrix< MT2, true > &rhs) |
Backend implementation of the Schur product between two column-major sparse matrices ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator% (const SparseMatrix< MT1, true > &lhs, const SparseMatrix< MT2, true > &rhs) |
Operator for the Schur product of two column-major sparse matrices ( ![]() | |
template<typename MT1 , typename MT2 > | |
decltype(auto) | blaze::operator- (const SparseMatrix< MT1, true > &lhs, const SparseMatrix< MT2, true > &rhs) |
Subtraction operator for the subtraction of two column-major sparse matrices ( ![]() | |
SparseMatrix global functions | |
template<typename MT , bool SO> | |
MT::Iterator | blaze::find (SparseMatrix< MT, SO > &sm, size_t i, size_t j) |
Searches for a specific matrix element. More... | |
template<typename MT , bool SO> | |
MT::ConstIterator | blaze::find (const SparseMatrix< MT, SO > &sm, size_t i, size_t j) |
Searches for a specific matrix element. More... | |
template<typename MT , bool SO> | |
MT::Iterator | blaze::lowerBound (SparseMatrix< MT, SO > &sm, size_t i, size_t j) |
Returns an iterator to the first index not less than the given index. More... | |
template<typename MT , bool SO> | |
MT::ConstIterator | blaze::lowerBound (const SparseMatrix< MT, SO > &sm, size_t i, size_t j) |
Returns an iterator to the first index not less than the given index. More... | |
template<typename MT , bool SO> | |
MT::Iterator | blaze::upperBound (SparseMatrix< MT, SO > &sm, size_t i, size_t j) |
Returns an iterator to the first index greater than the given index. More... | |
template<typename MT , bool SO> | |
MT::ConstIterator | blaze::upperBound (const SparseMatrix< MT, SO > &sm, size_t i, size_t j) |
Returns an iterator to the first index greater than the given index. More... | |
MatrixAccessProxy global functions | |
template<typename MT > | |
void | blaze::swap (const MatrixAccessProxy< MT > &a, const MatrixAccessProxy< MT > &b) noexcept |
Swapping the contents of two access proxies. More... | |
template<typename MT , typename T > | |
void | blaze::swap (const MatrixAccessProxy< MT > &a, T &b) noexcept |
Swapping the contents of an access proxy with another element. More... | |
template<typename T , typename MT > | |
void | blaze::swap (T &a, const MatrixAccessProxy< MT > &b) noexcept |
Swapping the contents of an access proxy with another element. More... | |
SparseMatrix operators | |
template<typename MT , bool SO, typename ST > | |
auto | blaze::operator*= (SparseMatrix< MT, SO > &mat, ST scalar) -> EnableIf_t< IsScalar_v< ST >, MT & > |
Multiplication assignment operator for the multiplication of a sparse matrix and a scalar value ( ![]() | |
template<typename MT , bool SO, typename ST > | |
auto | blaze::operator*= (SparseMatrix< MT, SO > &&mat, ST scalar) -> EnableIf_t< IsScalar_v< ST >, MT & > |
Multiplication assignment operator for the multiplication of a temporary sparse matrix and a scalar value ( ![]() | |
template<typename MT , bool SO, typename ST > | |
auto | blaze::operator/= (SparseMatrix< MT, SO > &mat, ST scalar) -> EnableIf_t< IsScalar_v< ST >, MT & > |
Division assignment operator for the division of a sparse matrix by a scalar value ( ![]() | |
template<typename MT , bool SO, typename ST > | |
auto | blaze::operator/= (SparseMatrix< MT, SO > &&mat, ST scalar) -> EnableIf_t< IsScalar_v< ST >, MT & > |
Division assignment operator for the division of a temporary sparse matrix by a scalar value ( ![]() | |
SparseMatrix functions | |
template<typename MT , bool SO> | |
bool | blaze::isnan (const SparseMatrix< MT, SO > &sm) |
Checks the given sparse matrix for not-a-number elements. More... | |
template<typename MT , bool SO> | |
bool | blaze::isinf (const SparseMatrix< MT, SO > &sm) |
Checks the given sparse matrix for infinite elements. More... | |
template<typename MT , bool SO> | |
bool | blaze::isfinite (const SparseMatrix< MT, SO > &sm) |
Checks the given sparse matrix for finite elements. More... | |
template<RelaxationFlag RF, typename MT , bool SO> | |
bool | blaze::isSymmetric (const SparseMatrix< MT, SO > &sm) |
Checks if the given sparse matrix is symmetric. More... | |
template<RelaxationFlag RF, typename MT , bool SO> | |
bool | blaze::isHermitian (const SparseMatrix< MT, SO > &sm) |
Checks if the given sparse matrix is Hermitian. More... | |
template<RelaxationFlag RF, typename MT , bool SO> | |
bool | blaze::isUniform (const SparseMatrix< MT, SO > &sm) |
Checks if the given sparse matrix is a uniform matrix. More... | |
template<RelaxationFlag RF, typename MT , bool SO> | |
bool | blaze::isZero (const SparseMatrix< MT, SO > &sm) |
Checks if the given sparse matrix is a zero matrix. More... | |
template<RelaxationFlag RF, typename MT , bool SO> | |
bool | blaze::isLower (const SparseMatrix< MT, SO > &sm) |
Checks if the given sparse matrix is a lower triangular matrix. More... | |
template<RelaxationFlag RF, typename MT , bool SO> | |
bool | blaze::isUniLower (const SparseMatrix< MT, SO > &sm) |
Checks if the given sparse matrix is a lower unitriangular matrix. More... | |
template<RelaxationFlag RF, typename MT , bool SO> | |
bool | blaze::isStrictlyLower (const SparseMatrix< MT, SO > &sm) |
Checks if the given sparse matrix is a strictly lower triangular matrix. More... | |
template<RelaxationFlag RF, typename MT , bool SO> | |
bool | blaze::isUpper (const SparseMatrix< MT, SO > &sm) |
Checks if the given sparse matrix is an upper triangular matrix. More... | |
template<RelaxationFlag RF, typename MT , bool SO> | |
bool | blaze::isUniUpper (const SparseMatrix< MT, SO > &sm) |
Checks if the given sparse matrix is an upper unitriangular matrix. More... | |
template<RelaxationFlag RF, typename MT , bool SO> | |
bool | blaze::isStrictlyUpper (const SparseMatrix< MT, SO > &sm) |
Checks if the given sparse matrix is a strictly upper triangular matrix. More... | |
template<RelaxationFlag RF, typename MT , bool SO> | |
bool | blaze::isDiagonal (const SparseMatrix< MT, SO > &sm) |
Checks if the give sparse matrix is diagonal. More... | |
template<RelaxationFlag RF, typename MT , bool SO> | |
bool | blaze::isIdentity (const SparseMatrix< MT, SO > &sm) |
Checks if the give sparse matrix is an identity matrix. More... | |
|
inline |
Applies the abs() function to each non-zero element of the sparse matrix sm.
sm | The input matrix. |
This function applies the abs() function to each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the abs() function:
|
inline |
Computes the inverse cosine for each non-zero element of the sparse matrix sm.
sm | The input matrix; all non-zero elements must be in the range ![]() |
The acos() function computes the inverse cosine for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the acos() function:
|
inline |
Computes the inverse hyperbolic cosine for each non-zero element of the sparse matrix sm.
sm | The input matrix; all non-zero elements must be in the range ![]() |
The acosh() function computes the inverse hyperbolic cosine for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the acosh() function:
|
inline |
Returns a matrix containing the phase angle of each single element of sm.
sm | The input matrix. |
The arg() function calculates the phase angle of each element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the arg() function:
|
inline |
Computes the inverse sine for each non-zero element of the sparse matrix sm.
sm | The input matrix; all non-zero elements must be in the range ![]() |
The asin() function computes the inverse sine for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the asin() function:
|
inline |
Computes the inverse hyperbolic sine for each non-zero element of the sparse matrix sm.
sm | The input matrix. |
The asinh() function computes the inverse hyperbolic sine for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the asinh() function:
|
inline |
Computes the inverse tangent for each non-zero element of the sparse matrix sm.
sm | The input matrix. |
The atan() function computes the inverse tangent for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the atan() function:
|
inline |
Computes the inverse hyperbolic tangent for each non-zero element of the sparse matrix sm.
sm | The input matrix; all non-zero elements must be in the range ![]() |
The atanh() function computes the inverse hyperbolic tangent for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the atanh() function:
|
inline |
Computes the cubic root of each non-zero element of the sparse matrix sm.
sm | The input matrix; all non-zero elements must be in the range ![]() |
The cbrt() function computes the cubic root of each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the cbrt() function:
|
inline |
Applies the ceil() function to each non-zero element of the sparse matrix sm.
sm | The input matrix. |
This function applies the ceil() function to each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the ceil() function:
|
inline |
Restricts each single element of the sparse matrix sm to the range .
sm | The input matrix. |
min | The lower delimiter. |
max | The upper delimiter. |
The clamp() function restricts each element of the input matrix sm to the range . The function returns an expression representing this operation.
The following example demonstrates the use of the clamp() function:
|
inline |
Returns a matrix containing the complex conjugate of each single element of sm.
sm | The input matrix. |
The conj() function calculates the complex conjugate of each element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the conj() function:
|
inline |
Computes the cosine for each non-zero element of the sparse matrix sm.
sm | The input matrix. |
The cos() function computes the cosine for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the cos() function:
|
inline |
Computes the hyperbolic cosine for each non-zero element of the sparse matrix sm.
sm | The input matrix. |
The cosh() function computes the hyperbolic cosine for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the cosh() function:
|
inline |
Returns the conjugate transpose matrix of sm.
sm | The input matrix. |
The ctrans() function returns an expression representing the conjugate transpose (also called adjoint matrix, Hermitian conjugate matrix or transjugate matrix) of the given input matrix sm.
The following example demonstrates the use of the ctrans() function:
Note that the ctrans() function has the same effect as manually applying the conj() and trans function in any order:
|
inline |
Declares the given sparse matrix expression sm as diagonal.
sm | The input matrix. |
std::invalid_argument | Invalid diagonal matrix specification. |
The decldiag function declares the given sparse matrix expression sm as diagonal. In case the given matrix is not a square matrix, a std::invalid_argument exception is thrown.
The following example demonstrates the use of the decldiag function:
|
inline |
Declares the given sparse matrix expression sm as Hermitian.
sm | The input matrix. |
std::invalid_argument | Invalid Hermitian matrix specification. |
The declherm function declares the given sparse matrix expression sm as Hermitian. In case the given matrix is not a square matrix, a std::invalid_argument exception is thrown.
The following example demonstrates the use of the declherm function:
|
inline |
Declares the given sparse matrix expression sm as lower.
sm | The input matrix. |
std::invalid_argument | Invalid lower matrix specification. |
The decllow function declares the given sparse matrix expression sm as lower. In case the given matrix is not a square matrix, a std::invalid_argument exception is thrown.
The following example demonstrates the use of the decllow function:
|
inline |
Declares the given sparse matrix expression sm as strictly lower.
sm | The input matrix. |
std::invalid_argument | Invalid strictly lower matrix specification. |
The declstrlow function declares the given sparse matrix expression sm as strictly lower. In case the given matrix is not a square matrix, a std::invalid_argument exception is thrown.
The following example demonstrates the use of the declstrlow function:
|
inline |
Declares the given sparse matrix expression sm as strictly upper.
sm | The input matrix. |
std::invalid_argument | Invalid strictly upper matrix specification. |
The declstrupp function declares the given sparse matrix expression sm as strictly upper. In case the given matrix is not a square matrix, a std::invalid_argument exception is thrown.
The following example demonstrates the use of the declstrupp function:
|
inline |
Declares the given sparse matrix expression sm as symmetric.
sm | The input matrix. |
std::invalid_argument | Invalid symmetric matrix specification. |
The declsym function declares the given sparse matrix expression sm as symmetric. In case the given matrix is not a square matrix, a std::invalid_argument exception is thrown.
The following example demonstrates the use of the declsym function:
|
inline |
Declares the given sparse matrix expression sm as unilower.
sm | The input matrix. |
std::invalid_argument | Invalid unilower matrix specification. |
The declunilow function declares the given sparse matrix expression sm as unilower. In case the given matrix is not a square matrix, a std::invalid_argument exception is thrown.
The following example demonstrates the use of the declunilow function:
|
inline |
Declares the given sparse matrix expression sm as uniupper.
sm | The input matrix. |
std::invalid_argument | Invalid uniupper matrix specification. |
The decluniupp function declares the given sparse matrix expression sm as uniupper. In case the given matrix is not a square matrix, a std::invalid_argument exception is thrown.
The following example demonstrates the use of the decluniupp function:
|
inline |
Declares the given sparse matrix expression sm as upper.
sm | The input matrix. |
std::invalid_argument | Invalid upper matrix specification. |
The declupp function declares the given sparse matrix expression sm as upper. In case the given matrix is not a square matrix, a std::invalid_argument exception is thrown.
The following example demonstrates the use of the declupp function:
|
inline |
Computes the error function for each non-zero element of the sparse matrix sm.
sm | The input matrix. |
The erf() function computes the error function for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the erf() function:
|
inline |
Computes the complementary error function for each non-zero element of the sparse matrix sm.
sm | The input matrix. |
The erfc() function computes the complementary error function for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the erfc() function:
|
inline |
Forces the evaluation of the given sparse matrix expression sm.
sm | The input matrix. |
The eval function forces the evaluation of the given sparse matrix expression sm. The function returns an expression representing the operation.
The following example demonstrates the use of the eval function:
|
inline |
Computes for each non-zero element of the sparse matrix sm.
sm | The input matrix. |
The exp() function computes for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the exp() function:
|
inline |
Computes for each non-zero element of the sparse matrix sm.
sm | The input matrix. |
The exp10() function computes for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the exp10() function:
|
inline |
Computes for each non-zero element of the sparse matrix sm.
sm | The input matrix. |
The exp2() function computes for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the exp2() function:
MT::ConstIterator blaze::find | ( | const SparseMatrix< MT, SO > & | sm, |
size_t | i, | ||
size_t | j | ||
) |
Searches for a specific matrix element.
sm | The given sparse matrix. |
i | The row index of the search element. The index has to be in the range ![]() |
j | The column index of the search element. The index has to be in the range ![]() |
This function can be used to check whether a specific element is contained in the sparse matrix. It specifically searches for the element with row index i and column index j. In case the element is found, the function returns an row/column iterator to the element. Otherwise an iterator just past the last non-zero element of row i or column j (the end() iterator) is returned. Note that the returned sparse matrix iterator is subject to invalidation due to inserting operations via the function call operator, the set() function or the insert() function!
MT::Iterator blaze::find | ( | SparseMatrix< MT, SO > & | sm, |
size_t | i, | ||
size_t | j | ||
) |
Searches for a specific matrix element.
sm | The given sparse matrix. |
i | The row index of the search element. The index has to be in the range ![]() |
j | The column index of the search element. The index has to be in the range ![]() |
This function can be used to check whether a specific element is contained in the sparse matrix. It specifically searches for the element with row index i and column index j. In case the element is found, the function returns an row/column iterator to the element. Otherwise an iterator just past the last non-zero element of row i or column j (the end() iterator) is returned. Note that the returned sparse matrix iterator is subject to invalidation due to inserting operations via the function call operator, the set() function or the insert() function!
|
noexcept |
Fixing the size of the given sparse matrix.
sm | The sparse matrix to be size-fixed. |
This function returns an expression representing the size-fixed given sparse matrix:
|
inline |
Applies the floor() function to each non-zero element of the sparse matrix sm.
sm | The input matrix. |
This function applies the floor() function to each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the floor() function:
|
inline |
Evaluates the given custom operation on each non-zero element of the sparse matrix sm.
sm | The input matrix. |
op | The custom operation. |
The forEach() function evaluates the given custom operation on each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the forEach() function:
|
inline |
Returns a matrix containing the imaginary parts of each single element of sm.
sm | The input matrix. |
The imag() function calculates the imaginary part of each element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the imag() function:
|
inline |
Computes the inverse cubic root of each non-zero element of the sparse matrix sm.
sm | The input matrix; all non-zero elements must be in the range ![]() |
The invcbrt() function computes the inverse cubic root of each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the invcbrt() function:
|
inline |
Computes the inverse square root of each non-zero element of the sparse matrix sm.
sm | The input matrix; all non-zero elements must be in the range ![]() |
The invsqrt() function computes the inverse square root of each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the invsqrt() function:
bool blaze::isDiagonal | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks if the give sparse matrix is diagonal.
sm | The sparse matrix to be checked. |
This function tests whether the matrix is diagonal, i.e. if the non-diagonal elements are default elements. In case of integral or floating point data types, a diagonal matrix has the form
\f[\left(\begin{array}{*{5}{c}} aa & 0 & 0 & \cdots & 0 \\ 0 & bb & 0 & \cdots & 0 \\ 0 & 0 & cc & \cdots & 0 \\ \vdots & \vdots & \vdots & \ddots & 0 \\ 0 & 0 & 0 & 0 & xx \\ \end{array}\right)\f]
The following example demonstrates the use of the function:
Optionally, it is possible to switch between strict semantics (blaze::strict) and relaxed semantics (blaze::relaxed):
It is also possible to check if a matrix expression results in a diagonal matrix:
However, note that this might require the complete evaluation of the expression, including the generation of a temporary matrix.
bool blaze::isfinite | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks the given sparse matrix for finite elements.
sm | The sparsre matrix to be checked for finite elements. |
This function checks if all elements of the sparse matrix are finite elements (i.e. normal, subnormal or zero elements, but not infinite or NaN). If all elements of the matrix are finite, the function returns true, otherwise it returns false.
bool blaze::isHermitian | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks if the given sparse matrix is Hermitian.
sm | The sparse matrix to be checked. |
This function checks if the given sparse matrix is an Hermitian matrix. The matrix is considered to be an Hermitian matrix if it is a square matrix whose conjugate transpose is equal to itself ( ), i.e. each matrix element
is equal to the complex conjugate of the element
. The following code example demonstrates the use of the function:
Optionally, it is possible to switch between strict semantics (blaze::strict) and relaxed semantics (blaze::relaxed):
It is also possible to check if a matrix expression results in an Hermitian matrix:
However, note that this might require the complete evaluation of the expression, including the generation of a temporary matrix.
bool blaze::isIdentity | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks if the give sparse matrix is an identity matrix.
sm | The sparse matrix to be checked. |
This function tests whether the matrix is an identity matrix, i.e. if the diagonal elements are 1 and the non-diagonal elements are 0. In case of integral or floating point data types, an identity matrix has the form
\f[\left(\begin{array}{*{5}{c}} 1 & 0 & 0 & \cdots & 0 \\ 0 & 1 & 0 & \cdots & 0 \\ 0 & 0 & 1 & \cdots & 0 \\ \vdots & \vdots & \vdots & \ddots & 0 \\ 0 & 0 & 0 & 0 & 1 \\ \end{array}\right)\f]
The following example demonstrates the use of the function:
Optionally, it is possible to switch between strict semantics (blaze::strict) and relaxed semantics (blaze::relaxed):
It is also possible to check if a matrix expression results in an identity matrix:
However, note that this might require the complete evaluation of the expression, including the generation of a temporary matrix.
bool blaze::isinf | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks the given sparse matrix for infinite elements.
sm | The sparse matrix to be checked for infinite elements. |
This function checks the sparse matrix for infinite (inf) elements. If at least one element of the matrix is infinite, the function returns true, otherwise it returns false.
bool blaze::isLower | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks if the given sparse matrix is a lower triangular matrix.
sm | The sparse matrix to be checked. |
This function checks if the given sparse matrix is a lower triangular matrix. The matrix is considered to be lower triangular if it is a square matrix of the form
\f[\left(\begin{array}{*{5}{c}} l_{0,0} & 0 & 0 & \cdots & 0 \\ l_{1,0} & l_{1,1} & 0 & \cdots & 0 \\ l_{2,0} & l_{2,1} & l_{2,2} & \cdots & 0 \\ \vdots & \vdots & \vdots & \ddots & \vdots \\ l_{N,0} & l_{N,1} & l_{N,2} & \cdots & l_{N,N} \\ \end{array}\right).\f]
or
matrices are considered as trivially lower triangular. The following code example demonstrates the use of the function:
Optionally, it is possible to switch between strict semantics (blaze::strict) and relaxed semantics (blaze::relaxed):
It is also possible to check if a matrix expression results in a lower triangular matrix:
However, note that this might require the complete evaluation of the expression, including the generation of a temporary matrix.
bool blaze::isnan | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks the given sparse matrix for not-a-number elements.
sm | The sparse matrix to be checked for not-a-number elements. |
This function checks the sparse matrix for not-a-number (NaN) elements. If at least one element of the matrix is not-a-number, the function returns true, otherwise it returns false.
bool blaze::isStrictlyLower | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks if the given sparse matrix is a strictly lower triangular matrix.
sm | The sparse matrix to be checked. |
This function checks if the given sparse matrix is a strictly lower triangular matrix. The matrix is considered to be strictly lower triangular if it is a square matrix of the form
\f[\left(\begin{array}{*{5}{c}} 0 & 0 & 0 & \cdots & 0 \\ l_{1,0} & 0 & 0 & \cdots & 0 \\ l_{2,0} & l_{2,1} & 0 & \cdots & 0 \\ \vdots & \vdots & \vdots & \ddots & \vdots \\ l_{N,0} & l_{N,1} & l_{N,2} & \cdots & 0 \\ \end{array}\right).\f]
The following code example demonstrates the use of the function:
Optionally, it is possible to switch between strict semantics (blaze::strict) and relaxed semantics (blaze::relaxed):
It is also possible to check if a matrix expression results in a strictly lower triangular matrix:
However, note that this might require the complete evaluation of the expression, including the generation of a temporary matrix.
bool blaze::isStrictlyUpper | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks if the given sparse matrix is a strictly upper triangular matrix.
sm | The sparse matrix to be checked. |
This function checks if the given sparse matrix is a strictly upper triangular matrix. The matrix is considered to be strictly upper triangular if it is a square matrix of the form
\f[\left(\begin{array}{*{5}{c}} 0 & u_{0,1} & u_{0,2} & \cdots & u_{0,N} \\ 0 & 0 & u_{1,2} & \cdots & u_{1,N} \\ 0 & 0 & 0 & \cdots & u_{2,N} \\ \vdots & \vdots & \vdots & \ddots & \vdots \\ 0 & 0 & 0 & \cdots & 0 \\ \end{array}\right).\f]
The following code example demonstrates the use of the function:
Optionally, it is possible to switch between strict semantics (blaze::strict) and relaxed semantics (blaze::relaxed):
It is also possible to check if a matrix expression results in a strictly upper triangular matrix:
However, note that this might require the complete evaluation of the expression, including the generation of a temporary matrix.
bool blaze::isSymmetric | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks if the given sparse matrix is symmetric.
sm | The sparse matrix to be checked. |
This function checks if the given sparse matrix is symmetric. The matrix is considered to be symmetric if it is a square matrix whose transpose is equal to itself ( ). The following code example demonstrates the use of the function:
Optionally, it is possible to switch between strict semantics (blaze::strict) and relaxed semantics (blaze::relaxed):
It is also possible to check if a matrix expression results in a symmetric matrix:
However, note that this might require the complete evaluation of the expression, including the generation of a temporary matrix.
bool blaze::isUniform | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks if the given sparse matrix is a uniform matrix.
sm | The sparse matrix to be checked. |
This function checks if the given sparse matrix is a uniform matrix. The matrix is considered to be uniform if all its elements are identical. The following code example demonstrates the use of the function:
Optionally, it is possible to switch between strict semantics (blaze::strict) and relaxed semantics (blaze::relaxed):
It is also possible to check if a matrix expression results in a uniform matrix:
However, note that this might require the complete evaluation of the expression, including the generation of a temporary matrix.
bool blaze::isUniLower | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks if the given sparse matrix is a lower unitriangular matrix.
sm | The sparse matrix to be checked. |
This function checks if the given sparse matrix is a lower unitriangular matrix. The matrix is considered to be lower unitriangular if it is a square matrix of the form
\f[\left(\begin{array}{*{5}{c}} 1 & 0 & 0 & \cdots & 0 \\ l_{1,0} & 1 & 0 & \cdots & 0 \\ l_{2,0} & l_{2,1} & 1 & \cdots & 0 \\ \vdots & \vdots & \vdots & \ddots & \vdots \\ l_{N,0} & l_{N,1} & l_{N,2} & \cdots & 1 \\ \end{array}\right).\f]
The following code example demonstrates the use of the function:
Optionally, it is possible to switch between strict semantics (blaze::strict) and relaxed semantics (blaze::relaxed):
It is also possible to check if a matrix expression results in a lower unitriangular matrix:
However, note that this might require the complete evaluation of the expression, including the generation of a temporary matrix.
bool blaze::isUniUpper | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks if the given sparse matrix is an upper unitriangular matrix.
sm | The sparse matrix to be checked. |
This function checks if the given sparse matrix is an upper unitriangular matrix. The matrix is considered to be upper unitriangular if it is a square matrix of the form
\f[\left(\begin{array}{*{5}{c}} 1 & u_{0,1} & u_{0,2} & \cdots & u_{0,N} \\ 0 & 1 & u_{1,2} & \cdots & u_{1,N} \\ 0 & 0 & 1 & \cdots & u_{2,N} \\ \vdots & \vdots & \vdots & \ddots & \vdots \\ 0 & 0 & 0 & \cdots & 1 \\ \end{array}\right).\f]
The following code example demonstrates the use of the function:
Optionally, it is possible to switch between strict semantics (blaze::strict) and relaxed semantics (blaze::relaxed):
It is also possible to check if a matrix expression results in an upper unitriangular matrix:
However, note that this might require the complete evaluation of the expression, including the generation of a temporary matrix.
bool blaze::isUpper | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks if the given sparse matrix is an upper triangular matrix.
sm | The sparse matrix to be checked. |
This function checks if the given sparse matrix is an upper triangular matrix. The matrix is considered to be upper triangular if it is a square matrix of the form
\f[\left(\begin{array}{*{5}{c}} u_{0,0} & u_{0,1} & u_{0,2} & \cdots & u_{0,N} \\ 0 & u_{1,1} & u_{1,2} & \cdots & u_{1,N} \\ 0 & 0 & u_{2,2} & \cdots & u_{2,N} \\ \vdots & \vdots & \vdots & \ddots & \vdots \\ 0 & 0 & 0 & \cdots & u_{N,N} \\ \end{array}\right).\f]
or
matrices are considered as trivially upper triangular. The following code example demonstrates the use of the function:
Optionally, it is possible to switch between strict semantics (blaze::strict) and relaxed semantics (blaze::relaxed):
It is also possible to check if a matrix expression results in an upper triangular matrix:
However, note that this might require the complete evaluation of the expression, including the generation of a temporary matrix.
bool blaze::isZero | ( | const SparseMatrix< MT, SO > & | sm | ) |
Checks if the given sparse matrix is a zero matrix.
sm | The sparse matrix to be checked. |
This function checks if the given sparse matrix is a zero matrix. The matrix is considered to be zero if all its elements are zero. The following code example demonstrates the use of the function:
Optionally, it is possible to switch between strict semantics (blaze::strict) and relaxed semantics (blaze::relaxed):
It is also possible to check if a matrix expression results is a zero matrix:
However, note that this might require the complete evaluation of the expression, including the generation of a temporary matrix.
|
inline |
Computes the Kronecker product of a dense matrix and a sparse matrix ( ).
lhs | The left-hand side dense matrix for the Kronecker product. |
rhs | The right-hand side sparse matrix for the Kronecker product. |
This kron() function computes the Kronecker product of the given dense matrix and sparse matrix:
The function returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
|
inline |
Computes the Kronecker product of two row-major sparse matrices ( ).
lhs | The left-hand side sparse matrix for the Kronecker product. |
rhs | The right-hand side sparse matrix for the Kronecker product. |
The kron() function computes the Kronecker product of the two given row-major sparse matrices:
The function returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
|
inline |
Operator for the Kronecker product of a row-major and a column-major sparse matrix ( ).
lhs | The left-hand side sparse matrix for the Kronecker product. |
rhs | The right-hand side sparse matrix for the Kronecker product. |
The kron() function computes the Kronecker product for a row-major sparse matrix and a column-major sparse matrix:
The function returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
|
inline |
Computes the Kronecker product of a sparse matrix and a dense matrix ( ).
lhs | The left-hand side sparse matrix for the Kronecker product. |
rhs | The right-hand side dense matrix for the Kronecker product. |
This kron() function computes the Kronecker product of the given sparse matrix and dense matrix:
The function returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
|
inline |
Operator for the Kronecker product of a column-major and a row-major sparse matrix ( ).
lhs | The left-hand side sparse matrix for the Kronecker product. |
rhs | The right-hand side sparse matrix for the Kronecker product. |
The kron() function computes the Kronecker product for a column-major sparse matrix and a row-major sparse matrix:
The function returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
|
inline |
Computes the Kronecker product of two column-major sparse matrices ( ).
lhs | The left-hand side sparse matrix for the Kronecker product. |
rhs | The right-hand side sparse matrix for the Kronecker product. |
The kron() function computes the Kronecker product of the two given column-major sparse matrices:
The function returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
decltype(auto) blaze::l1Norm | ( | const SparseMatrix< MT, SO > & | sm | ) |
Computes the L1 norm for the given sparse matrix.
sm | The given sparse matrix for the norm computation. |
This function computes the L1 norm of the given sparse matrix:
decltype(auto) blaze::l2Norm | ( | const SparseMatrix< MT, SO > & | sm | ) |
Computes the L2 norm for the given sparse matrix.
sm | The given sparse matrix for the norm computation. |
This function computes the L2 norm of the given sparse matrix:
decltype(auto) blaze::l3Norm | ( | const SparseMatrix< MT, SO > & | sm | ) |
Computes the L3 norm for the given sparse matrix.
sm | The given sparse matrix for the norm computation. |
This function computes the L3 norm of the given sparse matrix:
decltype(auto) blaze::l4Norm | ( | const SparseMatrix< MT, SO > & | sm | ) |
Computes the L4 norm for the given sparse matrix.
sm | The given sparse matrix for the norm computation. |
This function computes the L4 norm of the given sparse matrix:
|
inline |
Computes the natural logarithm of the absolute value of the gamma function for each non-zero element of the sparse matrix sm.
sm | The input matrix; all elements must be in the range ![]() |
The lgamma() function computes the natural logarithm of the absolute value of the gamma function for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the lgamma() function:
decltype(auto) blaze::linfNorm | ( | const SparseMatrix< MT, SO > & | sm | ) |
Computes the infinity norm for the given sparse matrix.
sm | The given sparse matrix for the norm computation. |
This function computes the infinity norm of the given sparse matrix:
|
inline |
Computes the natural logarithm for each non-zero element of the sparse matrix sm.
sm | The input matrix; all non-zero elements must be in the range ![]() |
The log() function computes the natural logarithm for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the log() function:
|
inline |
Computes the common logarithm for each non-zero element of the sparse matrix sm.
sm | The input matrix; all non-zero elements must be in the range ![]() |
The log10() function computes the common logarithm for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the log10() function:
|
inline |
Computes the natural logarithm of x+1 for each non-zero element of the sparse matrix sm.
sm | The input matrix; all elements must be in the range ![]() |
The log1p() function computes the natural logarithm of x+1 for each non-zero element of the input matrix sm. This may be preferred over the natural logarithm for higher precision computing the natural logarithm of a quantity very close to 1. The function returns an expression representing this operation.
The following example demonstrates the use of the log1p() function:
|
inline |
Computes the binary logarithm for each non-zero element of the sparse matrix sm.
sm | The input matrix; all non-zero elements must be in the range ![]() |
The log2() function computes the binary logarithm for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the log2() function:
MT::ConstIterator blaze::lowerBound | ( | const SparseMatrix< MT, SO > & | sm, |
size_t | i, | ||
size_t | j | ||
) |
Returns an iterator to the first index not less than the given index.
sm | The given sparse matrix. |
i | The row index of the search element. The index has to be in the range ![]() |
j | The column index of the search element. The index has to be in the range ![]() |
In case of a row-major matrix, this function returns a row iterator to the first element with an index not less than the given column index. In case of a column-major matrix, the function returns a column iterator to the first element with an index not less than the given row index. In combination with the upperBound() function this function can be used to create a pair of iterators specifying a range of indices. Note that the returned sparse matrix iterator is subject to invalidation due to inserting operations via the function call operator, the set() function or the insert() function!
MT::Iterator blaze::lowerBound | ( | SparseMatrix< MT, SO > & | sm, |
size_t | i, | ||
size_t | j | ||
) |
Returns an iterator to the first index not less than the given index.
sm | The given sparse matrix. |
i | The row index of the search element. The index has to be in the range ![]() |
j | The column index of the search element. The index has to be in the range ![]() |
In case of a row-major matrix, this function returns a row iterator to the first element with an index not less than the given column index. In case of a column-major matrix, the function returns a column iterator to the first element with an index not less than the given row index. In combination with the upperBound() function this function can be used to create a pair of iterators specifying a range of indices. Note that the returned sparse matrix iterator is subject to invalidation due to inserting operations via the function call operator, the set() function or the insert() function!
|
inline |
Computes the Lp norm for the given sparse matrix.
sm | The given sparse matrix for the norm computation. |
This function computes the Lp norm of the given sparse matrix, where the norm is specified by the runtime argument P:
decltype(auto) blaze::lpNorm | ( | const SparseMatrix< MT, SO > & | sm, |
ST | p | ||
) |
Computes the Lp norm for the given sparse matrix.
sm | The given sparse matrix for the norm computation. |
p | The norm parameter (p > 0). |
This function computes the Lp norm of the given sparse matrix, where the norm is specified by the runtime argument p:
|
inline |
Evaluates the given custom operation on each non-zero element of the sparse matrix sm.
sm | The input matrix. |
op | The custom operation. |
The map() function evaluates the given custom operation on each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the map() function:
|
inline |
Returns the largest element of the sparse matrix.
sm | The given sparse matrix. |
This function returns the largest non-zero element of the given sparse matrix. This function can only be used for element types that support the smaller-than relationship. In case the given matrix currently has either 0 rows or 0 columns, the returned value is the default value (e.g. 0 in case of fundamental data types).
|
inline |
Returns the largest element of each row/columns of the sparse matrix.
sm | The given sparse matrix. |
This function returns the largest element of each row/column of the given sparse matrix sm. In case the reduction flag RF is set to blaze::columnwise, a row vector containing the largest element of each column is returned. In case RF is set to blaze::rowwise, a column vector containing the largest element of each row is returned.
decltype(auto) blaze::maxNorm | ( | const SparseMatrix< MT, SO > & | sm | ) |
Computes the maximum norm for the given sparse matrix.
sm | The given sparse matrix for the norm computation. |
This function computes the maximum norm of the given sparse matrix:
|
inline |
Computes the (arithmetic) mean for the given sparse matrix.
sm | The given sparse matrix for the mean computation. |
std::invalid_argument | Invalid input matrix. |
This function computes the (arithmetic) mean for the given sparse matrix sm. Both the non-zero and zero elements of the sparse matrix are taken into account. Example:
In case the number of rows or columns of the given matrix is 0, a std::invalid_argument is thrown.
decltype(auto) blaze::mean | ( | const SparseMatrix< MT, SO > & | sm | ) |
Computes the row-/columnwise mean function for the given sparse matrix.
sm | The given sparse matrix for the mean computation. |
std::invalid_argument | Invalid input matrix. |
This function computes the row-/columnwise (arithmetic) mean for the given sparse matrix sm. In case RF is set to rowwise, the function returns a column vector containing the mean of each row of sm. In case RF is set to columnwise, the function returns a row vector containing the mean of each column of sm. Both the non-zero and zero elements of the sparse matrix are taken into account. Example:
In case RF is set to rowwise and the number of columns of the given matrix is 0 or in case RF is set to columnwise and the number of rows of the given matrix is 0, a std::invalid_argument is thrown.
|
inline |
Returns the smallest element of the sparse matrix.
sm | The given sparse matrix. |
This function returns the smallest non-zero element of the given sparse matrix. This function can only be used for element types that support the smaller-than relationship. In case the given matrix currently has either 0 rows or 0 columns, the returned value is the default value (e.g. 0 in case of fundamental data types).
|
inline |
Returns the smallest element of each row/columns of the sparse matrix.
sm | The given sparse matrix. |
This function returns the smallest non-zero element of each row/column of the given sparse matrix sm. In case the reduction flag RF is set to blaze::columnwise, a row vector containing the smallest element of each column is returned. In case RF is set to blaze::rowwise, a column vector containing the smallest element of each row is returned.
decltype(auto) blaze::minNorm | ( | const SparseMatrix< MT, SO > & | sm | ) |
Computes the minimum norm for the given sparse matrix.
sm | The given sparse matrix for the norm computation. |
This function computes the minimum norm of the given sparse matrix:
|
inline |
Forces the non-aliased evaluation of the given sparse matrix expression sm.
sm | The input matrix. |
The noalias function forces the non-aliased evaluation of the given sparse matrix expression sm. The function returns an expression representing the operation.
The following example demonstrates the use of the noalias function:
decltype(auto) blaze::norm | ( | const SparseMatrix< MT, SO > & | sm | ) |
Computes the L2 norm for the given sparse matrix.
sm | The given sparse matrix for the norm computation. |
This function computes the L2 norm of the given sparse matrix:
|
inline |
Disables the SIMD evaluation of the given sparse matrix expression sm.
sm | The input matrix. |
The nosimd function disables the SIMD evaluation of the given sparse matrix expression sm. The function returns an expression representing the operation.
The following example demonstrates the use of the nosimd function:
|
inline |
Inequality operator for the comparison of two sparse matrices.
lhs | The left-hand side sparse matrix for the comparison. |
rhs | The right-hand side sparse matrix for the comparison. |
|
inline |
Operator for the Schur product of a row-major sparse matrix and a row-major dense matrix ( ).
lhs | The left-hand side sparse matrix for the Schur product. |
rhs | The right-hand side dense matrix for the Schur product. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the Schur product of a row-major sparse matrix and a row-major dense matrix:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
In case the current sizes of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Operator for the Schur product of a row-major sparse matrix and a column-major dense matrix ( ).
lhs | The left-hand side sparse matrix for the Schur product. |
rhs | The right-hand side dense matrix for the Schur product. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the Schur product of a row-major sparse matrix and a column-major dense matrix:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
In case the current sizes of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Operator for the Schur product of two row-major sparse matrices ( ).
lhs | The left-hand side sparse matrix for the Schur product. |
rhs | The right-hand side sparse matrix for the Schur product. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the Schur product of two row-major sparse matrices:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Operator for the Schur product of a row-major and a column-major sparse matrix ( ).
lhs | The left-hand side sparse matrix for the Schur product. |
rhs | The right-hand side sparse matrix for the Schur product. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the Schur product of a row-major and a column-major sparse matrix:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Operator for the Schur product of a column-major sparse matrix and a row-major dense matrix ( ).
lhs | The left-hand side sparse matrix for the Schur product. |
rhs | The right-hand side dense matrix for the Schur product. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the Schur product of a column-major sparse matrix and a row-major dense matrix:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
In case the current sizes of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Operator for the Schur product of a column-major sparse matrix and a column-major dense matrix ( ).
lhs | The left-hand side sparse matrix for the Schur product. |
rhs | The right-hand side dense matrix for the Schur product. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the Schur product of a column-major sparse matrix and a column-major dense matrix:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
In case the current sizes of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Operator for the Schur product of a column-major and a row-major sparse matrix ( ).
lhs | The left-hand side sparse matrix for the Schur product. |
rhs | The right-hand side sparse matrix for the Schur product. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the Schur product of a column-major and a row-major sparse matrix:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Operator for the Schur product of two column-major sparse matrices ( ).
lhs | The left-hand side sparse matrix for the Schur product. |
rhs | The right-hand side sparse matrix for the Schur product. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the Schur product of two column-major sparse matrices:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Multiplication operator for the dense vector-sparse vector outer product ( ).
lhs | The left-hand side dense vector for the outer product. |
rhs | The right-hand side transpose sparse vector for the outer product. |
This operator represents the outer product between a dense vector and a transpose sparse vector:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved element types VT1::ElementType and VT2::ElementType. Both vector types VT1 and VT2 as well as the two element types VT1::ElementType and VT2::ElementType have to be supported by the MultTrait class template.
|
inline |
Multiplication operator for the multiplication of a sparse matrix and a scalar value ( ).
mat | The left-hand side sparse matrix for the multiplication. |
scalar | The right-hand side scalar value for the multiplication. |
This operator represents the multiplication between a sparse matrix and a scalar value:
The operator returns an expression representing a sparse matrix of the higher-order element type of the involved data types MT::ElementType and T2. Note that this operator only works for scalar values of built-in data type.
|
inline |
Multiplication operator for the multiplication of two row-major sparse matrices ( ).
lhs | The left-hand side sparse matrix for the matrix multiplication. |
rhs | The right-hand side sparse matrix for the matrix multiplication. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the multiplication of two row-major sparse matrices:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Multiplication operator for the multiplication of a row-major sparse matrix and a column-major sparse matrix ( ).
lhs | The left-hand side sparse matrix for the matrix multiplication. |
rhs | The right-hand side sparse matrix for the matrix multiplication. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the multiplication of a row-major sparse matrix and a column-major sparse matrix:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Multiplication operator for the multiplication of a column-major sparse matrix and a row-major sparse matrix ( ).
lhs | The left-hand side sparse matrix for the matrix multiplication. |
rhs | The right-hand side sparse matrix for the matrix multiplication. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the multiplication of a column-major sparse matrix and a row-major sparse matrix:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Multiplication operator for the multiplication of two column-major sparse matrices ( ).
lhs | The left-hand side sparse matrix for the matrix multiplication. |
rhs | The right-hand side sparse matrix for the matrix multiplication. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the multiplication of two column-major sparse matrices:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the MultTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Multiplication operator for the sparse vector-dense vector outer product ( ).
lhs | The left-hand side sparse vector for the outer product. |
rhs | The right-hand side transpose dense vector for the outer product. |
This operator represents the outer product between a sparse vector and a transpose dense vector:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved element types VT1::ElementType and VT2::ElementType. Both vector types VT1 and VT2 as well as the two element types VT1::ElementType and VT2::ElementType have to be supported by the MultTrait class template.
|
inline |
Multiplication operator for the sparse vector-sparse vector outer product ( ).
lhs | The left-hand side sparse vector for the outer product. |
rhs | The right-hand side transpose sparse vector for the outer product. |
This operator represents the outer product between a sparse vector and a transpose sparse vector:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved element types VT1::ElementType and VT2::ElementType. Both vector types VT1 and VT2 as well as the two element types VT1::ElementType and VT2::ElementType have to be supported by the MultTrait class template.
|
inline |
Multiplication operator for the multiplication of a scalar value and a sparse matrix ( ).
scalar | The left-hand side scalar value for the multiplication. |
mat | The right-hand side sparse matrix for the multiplication. |
This operator represents the multiplication between a scalar value and a sparse matrix:
The operator returns an expression representing a sparse matrix of the higher-order element type of the involved data types ST and MT::ElementType. Note that this operator only works for scalar values of built-in data type.
|
inline |
Multiplication assignment operator for the multiplication of a temporary sparse matrix and a scalar value ( ).
mat | The left-hand side temporary sparse matrix for the multiplication. |
scalar | The right-hand side scalar value for the multiplication. |
std::invalid_argument | Invalid scaling of restricted matrix. |
In case the matrix MT is restricted and the assignment would violate an invariant of the matrix, a std::invalid_argument exception is thrown.
|
inline |
Multiplication assignment operator for the multiplication of a sparse matrix and a scalar value ( ).
mat | The left-hand side sparse matrix for the multiplication. |
scalar | The right-hand side scalar value for the multiplication. |
std::invalid_argument | Invalid scaling of restricted matrix. |
In case the matrix MT is restricted and the assignment would violate an invariant of the matrix, a std::invalid_argument exception is thrown.
|
inline |
Addition operator for the addition of two row-major sparse matrices ( ).
lhs | The left-hand side sparse matrix for the matrix addition. |
rhs | The right-hand side sparse matrix to be added to the left-hand side matrix. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the addition of two row-major sparse matrices:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the AddTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Addition operator for the addition of a row-major and a column-major sparse matrix ( ).
lhs | The left-hand side sparse matrix for the matrix addition. |
rhs | The right-hand side sparse matrix to be added to the left-hand side matrix. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the addition of a row-major and a column-major sparse matrix:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the AddTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Addition operator for the addition of a column-major and a row-major sparse matrix ( ).
lhs | The left-hand side sparse matrix for the matrix addition. |
rhs | The right-hand side sparse matrix to be added to the left-hand side matrix. |
std::invalid_argument | Matrix sizes do not match |
This operator represents the addition of a column-major and a row-major sparse matrix:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the AddTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Addition operator for the addition of two column-major sparse matrices ( ).
lhs | The left-hand side sparse matrix for the matrix addition. |
rhs | The right-hand side sparse matrix to be added to the left-hand side matrix. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the addition of two column-major sparse matrices:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the AddTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Unary minus operator for the negation of a sparse matrix ( ).
sm | The sparse matrix to be negated. |
This operator represents the negation of a sparse matrix:
The operator returns an expression representing the negation of the given sparse matrix.
|
inline |
Subtraction operator for the subtraction of two row-major sparse matrices ( ).
lhs | The left-hand side sparse matrix for the matrix subtraction. |
rhs | The right-hand side sparse matrix to be subtracted from the left-hand side matrix. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the subtraction of two row-major sparse matrices:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the SubTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Subtraction operator for the subtraction of a row-major and a column-major sparse matrix ( ).
lhs | The left-hand side sparse matrix for the matrix subtraction. |
rhs | The right-hand side sparse matrix to be subtracted from the left-hand side matrix. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the subtraction of a row-major and a column-major sparse matrix:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the SubTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Subtraction operator for the subtraction of a column-major and a row-major sparse matrix ( ).
lhs | The left-hand side sparse matrix for the matrix subtraction. |
rhs | The right-hand side sparse matrix to be subtracted from the left-hand side matrix. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the subtraction of a column-major and a row-major sparse matrix:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the SubTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Subtraction operator for the subtraction of two column-major sparse matrices ( ).
lhs | The left-hand side sparse matrix for the matrix subtraction. |
rhs | The right-hand side sparse matrix to be subtracted from the left-hand side matrix. |
std::invalid_argument | Matrix sizes do not match. |
This operator represents the subtraction of two column-major sparse matrices:
The operator returns an expression representing a sparse matrix of the higher-order element type of the two involved matrix element types MT1::ElementType and MT2::ElementType. Both matrix types MT1 and MT2 as well as the two element types MT1::ElementType and MT2::ElementType have to be supported by the SubTrait class template.
In case the current number of rows and columns of the two given matrices don't match, a std::invalid_argument is thrown.
|
inline |
Division operator for the division of a sparse matrix by a scalar value ( ).
mat | The left-hand side sparse matrix for the division. |
scalar | The right-hand side scalar value for the division. |
This operator represents the division of a sparse matrix by a scalar value:
The operator returns an expression representing a sparse matrix of the higher-order element type of the involved data types MT::ElementType and ST. Note that this operator only works for scalar values of built-in data type.
|
inline |
Division assignment operator for the division of a temporary sparse matrix by a scalar value ( ).
mat | The left-hand side temporary sparse matrix for the division. |
scalar | The right-hand side scalar value for the division. |
std::invalid_argument | Invalid scaling of restricted matrix. |
In case the matrix MT is restricted and the assignment would violate an invariant of the matrix, a std::invalid_argument exception is thrown.
|
inline |
Division assignment operator for the division of a sparse matrix by a scalar value ( ).
mat | The left-hand side sparse matrix for the division. |
scalar | The right-hand side scalar value for the division. |
std::invalid_argument | Invalid scaling of restricted matrix. |
In case the matrix MT is restricted and the assignment would violate an invariant of the matrix, a std::invalid_argument exception is thrown.
|
inline |
Equality operator for the comparison of two sparse matrices.
lhs | The left-hand side sparse matrix for the comparison. |
rhs | The right-hand side sparse matrix for the comparison. |
|
inline |
Computes the exponential value for each non-zero element of the sparse matrix sm.
sm | The input matrix. |
exp | The scalar exponent. |
The pow() function computes the exponential value for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the pow() function:
|
inline |
Reduces the given sparse matrix by means of multiplication.
sm | The given sparse matrix for the reduction operation. |
This function reduces the non-zero elements of the given sparse matrix sm by means of multiplication:
Please note that the evaluation order of the reduction operation is unspecified.
|
inline |
Reduces the given sparse matrix by means of multiplication.
sm | The given sparse matrix for the reduction operation. |
This function reduces the non-zero elements of the rows or columns of the given sparse matrix sm by means of multiplication. In case the reduction flag RF is set to blaze::columnwise, the elements of the matrix are reduced column-wise and the result is a row vector. In case RF is set to blaze::rowwise, the elements of the matrix are reduced row-wise and the result is a column vector:
Please note that the evaluation order of the reduction operation is unspecified.
|
inline |
Returns a matrix containing the real parts of each single element of sm.
sm | The input matrix. |
The real() function calculates the real part of each element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the real() function:
|
inline |
Performs a custom reduction operation on the given sparse matrix.
sm | The given sparse matrix for the reduction computation. |
op | The reduction operation. |
This function reduces the non-zero elements of the given sparse matrix sm by means of the given reduction operation op:
As demonstrated in the example it is possible to pass any binary callable as custom reduction operation. See Custom Operations for a detailed overview of the possibilities of custom operations.
Please note that the evaluation order of the reduction operation is unspecified. Thus the behavior is non-deterministic if op is not associative or not commutative. Also, the operation is undefined if the given reduction operation modifies the values.
|
inline |
Performs a custom reduction operation on the given sparse matrix.
sm | The given sparse matrix for the reduction computation. |
op | The reduction operation. |
This function reduces the rows or columns of the given sparse matrix sm by means of the given reduction operation op. In case the reduction flag RF is set to blaze::columnwise, the elements of the matrix are reduced column-wise and the result is a row vector. In case RF is set to blaze::rowwise, the elements of the matrix are reduced row-wise and the result is a column vector:
Please note that the evaluation order of the reduction operation is unspecified. Thus the behavior is non-deterministic if op is not associative or not commutative. Also, the operation is undefined if the given reduction operation modifies the values.
|
inline |
Repeats the given sparse matrix.
sm | The sparse matrix to be repeated. |
This function returns an expression representing the repeated sparse matrix:
|
inline |
Repeats the given sparse matrix.
sm | The sparse matrix to be repeated. |
m | The number of row-wise repetitions. |
n | The number of column-wise repetitions. |
This function returns an expression representing the repeated sparse matrix:
|
inline |
Applies the round() function to each non-zero element of the sparse matrix sm.
sm | The input matrix. |
This function applies the round() function to each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the round() function:
|
inline |
Forces the serial evaluation of the given sparse matrix expression sm.
sm | The input matrix. |
The serial function forces the serial evaluation of the given sparse matrix expression sm. The function returns an expression representing the operation.
The following example demonstrates the use of the serial function:
|
inline |
Applies the sign() function to each non-zero element of the sparse matrix sm.
sm | The input matrix. |
This function applies the sign() function to each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the sign() function:
|
inline |
Computes the sine for each non-zero element of the sparse matrix sm.
sm | The input matrix. |
The sin() function computes the sine for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the sin() function:
|
inline |
Computes the hyperbolic sine for each non-zero element of the sparse matrix sm.
sm | The input matrix. |
The sinh() function computes the hyperbolic sine for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the sinh() function:
decltype(auto) blaze::sqrNorm | ( | const SparseMatrix< MT, SO > & | sm | ) |
Computes the squared L2 norm for the given sparse matrix.
sm | The given sparse matrix for the norm computation. |
This function computes the squared L2 norm of the given sparse matrix:
|
inline |
Computes the square root of each non-zero element of the sparse matrix sm.
sm | The input matrix; all non-zero elements must be in the range ![]() |
The sqrt() function computes the square root of each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the sqrt() function:
|
inline |
Computes the standard deviation for the given sparse matrix.
sm | The given sparse matrix for the standard deviation computation. |
std::invalid_argument | Invalid input matrix. |
This function computes the standard deviation for the given sparse matrix sm. Both the non-zero and zero elements of the sparse matrix are taken into account. Example:
In case the size of the given matrix is smaller than 2, a std::invalid_argument is thrown.
|
inline |
Computes the row-/columnwise standard deviation function for the given sparse matrix.
sm | The given sparse matrix for the standard deviation computation. |
std::invalid_argument | Invalid input matrix. |
This function computes the row-/columnwise standard deviation for the given sparse matrix sm. In case RF is set to rowwise, the function returns a column vector containing the standard deviation of each row of sm. In case RF is set to columnwise, the function returns a row vector containing the standard deviation of each column of sm. Both the non-zero and zero elements of the sparse matrix are taken into account. Example:
In case RF is set to rowwise and the number of columns of the given matrix is smaller than 2 or in case RF is set to columnwise and the number of rows of the given matrix is smaller than 2, a std::invalid_argument is thrown.
|
inline |
Reduces the given sparse matrix by means of addition.
sm | The given sparse matrix for the reduction operation. |
This function reduces the non-zero elements of the given sparse matrix sm by means of addition:
Please note that the evaluation order of the reduction operation is unspecified.
|
inline |
Reduces the given sparse matrix by means of addition.
sm | The given sparse matrix for the reduction operation. |
This function reduces the non-zero elements of the rows or columns of the given sparse matrix sm by means of addition. In case the reduction flag RF is set to blaze::columnwise, the elements of the matrix are reduced column-wise and the result is a row vector. In case RF is set to blaze::rowwise, the elements of the matrix are reduced row-wise and the result is a column vector:
Please note that the evaluation order of the reduction operation is unspecified.
|
inlinenoexcept |
Swapping the contents of two access proxies.
a | The first access proxy to be swapped. |
b | The second access proxy to be swapped. |
|
inlinenoexcept |
Swapping the contents of an access proxy with another element.
a | The access proxy to be swapped. |
b | The other element to be swapped. |
|
inlinenoexcept |
Swapping the contents of an access proxy with another element.
a | The other element to be swapped. |
b | The access proxy to be swapped. |
|
inline |
Computes the tangent for each non-zero element of the sparse matrix sm.
sm | The input matrix. |
The tan() function computes the tangent for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the tan() function:
|
inline |
Computes the hyperbolic tangent for each non-zero element of the sparse matrix sm.
sm | The input matrix; all non-zero elements must be in the range ![]() |
The tanh() function computes the hyperbolic tangent for each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the tanh() function:
|
inline |
Calculation of the transpose of the given sparse matrix.
sm | The sparse matrix to be transposed. |
This function returns an expression representing the transpose of the given sparse matrix:
|
inline |
Conditional calculation of the transpose of the given sparse matrix.
sm | The sparse matrix to be transposed. |
In case the given compile time condition evaluates to true, this function returns an expression representing the transpose of the given sparse matrix. Otherwise the function returns a reference to the given matrix.
|
inline |
Applies the trunc() function to each non-zero element of the sparse matrix sm.
sm | The input matrix. |
This function applies the trunc() function to each non-zero element of the input matrix sm. The function returns an expression representing this operation.
The following example demonstrates the use of the trunc() function:
|
inline |
Backend implementation of the Schur product between two column-major sparse matrices ( ).
lhs | The left-hand side sparse matrix for the Schur product. |
rhs | The right-hand side sparse matrix for the Schur product. |
This function implements a performance optimized treatment of the Schur product between two column-major sparse matrices.
MT::ConstIterator blaze::upperBound | ( | const SparseMatrix< MT, SO > & | sm, |
size_t | i, | ||
size_t | j | ||
) |
Returns an iterator to the first index greater than the given index.
sm | The given sparse matrix. |
i | The row index of the search element. The index has to be in the range ![]() |
j | The column index of the search element. The index has to be in the range ![]() |
In case of a row-major matrix, this function returns a row iterator to the first element with an index greater than the given column index. In case of a column-major matrix, the function returns a column iterator to the first element with an index greater than the given row index. In combination with the lowerBound() function this function can be used to create a pair of iterators specifying a range of indices. Note that the returned sparse matrix iterator is subject to invalidation due to inserting operations via the function call operator, the set() function or or the insert() function!
MT::Iterator blaze::upperBound | ( | SparseMatrix< MT, SO > & | sm, |
size_t | i, | ||
size_t | j | ||
) |
Returns an iterator to the first index greater than the given index.
sm | The given sparse matrix. |
i | The row index of the search element. The index has to be in the range ![]() |
j | The column index of the search element. The index has to be in the range ![]() |
In case of a row-major matrix, this function returns a row iterator to the first element with an index greater than the given column index. In case of a column-major matrix, the function returns a column iterator to the first element with an index greater than the given row index. In combination with the lowerBound() function this function can be used to create a pair of iterators specifying a range of indices. Note that the returned sparse matrix iterator is subject to invalidation due to inserting operations via the function call operator, the set() function or or the insert() function!
|
inline |
Computes the variance for the given sparse matrix.
sm | The given sparse matrix for the variance computation. |
std::invalid_argument | Invalid input matrix. |
This function computes the variance for the given sparse matrix sm. Both the non-zero and zero elements of the sparse matrix are taken into account. Example:
In case the size of the given matrix is smaller than 2, a std::invalid_argument is thrown.
|
inline |
Computes the row-/column-wise variance function for the given sparse matrix.
sm | The given sparse matrix for the variance computation. |
std::invalid_argument | Invalid input matrix. |
This function computes the row-/column-wise variance for the given sparse matrix sm. In case RF is set to rowwise, the function returns a column vector containing the variance of each row of sm. In case RF is set to columnwise, the function returns a row vector containing the variance of each column of dm. Both the non-zero and zero elements of the sparse matrix are taken into account. Example:
In case RF is set to rowwise and the number of columns of the given matrix is smaller than 2 or in case RF is set to columnwise and the number of rows of the given matrix is smaller than 2, a std::invalid_argument is thrown.