Man page - rsbmatrix(3)

Packages contains this manual

Manual

RsbMatrix< NT >

NAME
SYNOPSIS
Public Types
Public Member Functions
Detailed Description
template<RSBP_Scalar_t NT>
Member Enumeration Documentation
template<RSBP_Scalar_t NT> enum RsbMatrix::RsbSym
Constructor & Destructor Documentation
template<RSBP_Scalar_t NT> RsbMatrix< NT >::RsbMatrix (rsb_coo_idx_t nrA,rsb_coo_idx_t ncA, const RsbSym sym = IsGen) [inline]
template<RSBP_Scalar_t NT> RsbMatrix< NT >::RsbMatrix (rsb_coo_idx_t nrA,const rsb_coo_idx_t * RP, const rsb_coo_idx_t * JA, const NT * VA,const RsbSym sym = IsGen) [inline]
template<RSBP_Scalar_t NT> RsbMatrix< NT >::RsbMatrix (const rsb_coo_idx_t* IA, const rsb_coo_idx_t * JA, const NT * VA, rsb_nnz_idx_t nnzA,const rsb_flags_t flagsA = RSB_FLAG_NOFLAGS) [inline]
template<RSBP_Scalar_t NT> RsbMatrix< NT >::RsbMatrix (const rsb_char_t *filename, const RsbSym sym = IsGen) [inline]
template<RSBP_Scalar_t NT> RsbMatrix< NT >::RsbMatrix (const RsbMatrix< NT> & A_Rsb, bool do_trans = false, rsb_flags_t flagsA =RSB_FLAG_NOFLAGS) [inline]
template<RSBP_Scalar_t NT> RsbMatrix< NT >::RsbMatrix (RsbMatrix< NT > &&other) [inline]
template<RSBP_Scalar_t NT> RsbMatrix< NT >::˜RsbMatrix (void) [inline]
Member Function Documentation
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::_add (rsb_coo_idx_t i, rsb_coo_idx_t j, NT val) [inline]
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::_close (void) [inline]
template<RSBP_Scalar_t NT> size_t RsbMatrix< NT >::_get_index_storage_bytes(void) const [inline]
template<RSBP_Scalar_t NT> size_t RsbMatrix< NT >::_get_storage_bytes(void) const [inline]
template<RSBP_Scalar_t NT> rsb_string_t RsbMatrix< NT >::_info (void) const[inline]
template<RSBP_Scalar_t NT> bool RsbMatrix< NT >::_is_complex (void) const[inline]
template<RSBP_Scalar_t NT> rsb_blk_idx_t RsbMatrix< NT >::blocks (void)const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::close (void)[inline]
template<RSBP_Scalar_t NT> rsb_coo_idx_t RsbMatrix< NT >::cols (void) const[inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::file_save (constrsb_char_t * filename = RSBP_NULL) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::get_coo (rsb_trans_t transA, NT * VA, rsb_coo_idx_t * IA,rsb_coo_idx_t * JA, rsb_flags_t flags) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::get_csr (rsb_trans_t transA, NT * VA, rsb_coo_idx_t * RP,rsb_coo_idx_t * JA, rsb_flags_t flags) const [inline]
template<RSBP_Scalar_t NT> rsb_flags_t RsbMatrix< NT >::get_flags_t (enumrsb_mif_t mif) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::get_info (enumrsb_mif_t miflags, void * minfop) const [inline]
template<RSBP_Scalar_t NT> rsb_blk_idx_t RsbMatrix< NT >::get_info_blk_t(enum rsb_mif_t mif) const [inline]
template<RSBP_Scalar_t NT> rsb_coo_idx_t RsbMatrix< NT >::get_info_coo_t(enum rsb_mif_t mif) const [inline]
template<RSBP_Scalar_t NT> rsb_nnz_idx_t RsbMatrix< NT >::get_info_nnz_t(enum rsb_mif_t mif) const [inline]
template<RSBP_Scalar_t NT> rsb_flags_t RsbMatrix< NT>::get_info_rsb_flags_t (enum rsb_mif_t mif) const [inline]
template<RSBP_Scalar_t NT> size_t RsbMatrix< NT >::get_info_size_t (enumrsb_mif_t mif) const [inline]
template<RSBP_Scalar_t NT> rsb_string_t RsbMatrix< NT >::get_info_str(const char * key) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::get_nrm (NT *Np, enum rsb_extff_t flags) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::get_rows_sparse (rsb_trans_t transA, const NT * alphap, NT * VA,rsb_coo_idx_t * IA, rsb_coo_idx_t * JA, rsb_coo_idx_t frA,rsb_coo_idx_t lrA, rsb_nnz_idx_t * rnzp, rsb_flags_t flags) const[inline]
template<RSBP_Scalar_t NT> rsb_type_t RsbMatrix< NT >::get_type_t (enumrsb_mif_t mif) const [inline]
template<RSBP_Scalar_t NT> NT RsbMatrix< NT >::get_val (const rsb_coo_idx_ti, const rsb_coo_idx_t j, rsb_flags_t flags = RSB_FLAG_NOFLAGS) const[inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::get_vals (NT *VA, const rsb_coo_idx_t * IA, const rsb_coo_idx_t * JA, rsb_nnz_idx_tnnz, rsb_flags_t flags) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::get_vec (NT *Dp, enum rsb_extff_t flags) const [inline]
template<RSBP_Scalar_t NT> rsb_nnz_idx_t RsbMatrix< NT >::nnz (void) const[inline]
template<RSBP_Scalar_t NT> NT RsbMatrix< NT >::normInf (void) const[inline]
template<RSBP_Scalar_t NT> NT RsbMatrix< NT >::normOne (void) const[inline]
template<RSBP_Scalar_t NT> bool RsbMatrix< NT >::operator!= (constRsbMatrix< NT > & B_Rsb) const [inline]
template<RSBP_Scalar_t NT> RsbMatrix & RsbMatrix< NT >::operator= (constRsbMatrix< NT > & A_Rsb) [inline]
template<RSBP_Scalar_t NT> bool RsbMatrix< NT >::operator== (constRsbMatrix< NT > & B_Rsb) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::rndr (constrsb_char_t * filename = RSBP_NULL, rsb_coo_idx_t pmWidth = 512,rsb_coo_idx_t pmHeight = 512, rsb_marf_t rflags = RSB_MARF_EPS) const[inline]
template<RSBP_Scalar_t NT> rsb_coo_idx_t RsbMatrix< NT >::rows (void) const[inline]
template<RSBP_Scalar_t NT> rsb_flags_t RsbMatrix< NT >::rsbflags (void)const [inline]
template<RSBP_Scalar_t NT> rsb_type_t RsbMatrix< NT >::rsbtype (void) const[inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::set_val (constNT val, const rsb_coo_idx_t i, const rsb_coo_idx_t j, rsb_flags_t flags= RSB_FLAG_NOFLAGS) [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::set_vals (constNT * VA, const rsb_coo_idx_t * IA, const rsb_coo_idx_t * JA,rsb_nnz_idx_t nnz, rsb_flags_t flags) [inline]
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::spmm (rsb_trans_t transA, const NT * alphap, rsb_coo_idx_t nrhs,rsb_flags_t order, const NT * Bp, rsb_nnz_idx_t ldB, const NT * betap,NT * Cp, rsb_nnz_idx_t ldC) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spmm(rsb_trans_t transA, const NT alpha, rsb_coo_idx_t nrhs, rsb_flags_torder, const NT * Bp, const NT beta, NT * Cp) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spmm(rsb_trans_t transA, const NT alpha, rsb_coo_idx_t nrhs, rsb_flags_torder, const NT * Bp, rsb_nnz_idx_t ldB, const NT beta, NT * Cp,rsb_nnz_idx_t ldC) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spmv (NT * y,const NT * x, bool do_trans = false) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::spmv (rsb_trans_t transA, const NT * alphap, const NT * Xp,rsb_coo_idx_t incX, const NT * betap, NT * Yp, rsb_coo_idx_t incY)const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spmv(rsb_trans_t transA, const NT alpha, const NT * Xp, const NT beta, NT *Yp) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spmv(rsb_trans_t transA, const NT alpha, const NT * Xp, rsb_coo_idx_t incX,const NT beta, NT * Yp, rsb_coo_idx_t incY) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsm (NT * y,const NT * x, rsb_coo_idx_t nrhs, bool do_trans = false) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsm (NT * y,rsb_coo_idx_t nrhs, bool do_trans = false) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::spsm (rsb_trans_t transT, const NT * alphap, rsb_coo_idx_t nrhs,rsb_flags_t order, const NT * betap, const NT * Bp, rsb_nnz_idx_t ldB,NT * Cp, rsb_nnz_idx_t ldC) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsm(rsb_trans_t transT, const NT alpha, rsb_coo_idx_t nrhs, const NT * Bp,NT * Cp) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsm(rsb_trans_t transT, const NT alpha, rsb_coo_idx_t nrhs, rsb_flags_torder, const NT beta, const NT * Bp, rsb_nnz_idx_t ldB, NT * Cp,rsb_nnz_idx_t ldC) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsv (NT * y,bool do_trans = false) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsv (NT * y,const NT * x, bool do_trans = false) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::spsv (rsb_trans_t transT, const NT * alphap, const NT * Xp,rsb_coo_idx_t incX, NT * Yp, rsb_coo_idx_t incY) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsv(rsb_trans_t transT, const NT alpha, const NT * Xp, NT * Yp) const[inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::tune_spmm(rsb_real_t & sf, rsb_trans_t transA, const NT alpha, rsb_coo_idx_tnrhs, rsb_flags_t order, const NT * Bp, const NT beta, NT * Cp)[inline]
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::tune_spmm (rsb_real_t * sfp, rsb_int_t * tnp, rsb_int_t maxr,rsb_time_t maxt, rsb_trans_t transA, const NT * alphap, rsb_coo_idx_tnrhs, rsb_flags_t order, const NT * Bp, rsb_nnz_idx_t ldB, const NT *betap, NT * Cp, rsb_nnz_idx_t ldC) [inline]
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::tune_spmm (rsb_real_t * sfp, rsb_int_t * tnp, rsb_int_t maxr,rsb_time_t maxt, rsb_trans_t transA, const NT alpha, rsb_coo_idx_tnrhs, rsb_flags_t order, const NT * Bp, rsb_nnz_idx_t ldB, const NTbeta, NT * Cp, rsb_nnz_idx_t ldC) [inline]
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::tune_spmm_threads (rsb_real_t * sfp = RSBP_NULL, rsb_int_t * tnp =RSBP_NULL, rsb_int_t maxr = 0, rsb_time_t maxt = 0, rsb_trans_t transA= RSB_TRANSPOSITION_N, const NT * alphap = RSBP_NULL, rsb_coo_idx_tnrhs = 1, rsb_flags_t order = RSB_FLAG_WANT_COLUMN_MAJOR_ORDER, constNT * Bp = RSBP_NULL, rsb_nnz_idx_t ldB = 0, const NT * betap =RSBP_NULL, NT * Cp = RSBP_NULL, rsb_nnz_idx_t ldC = 0) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::tune_spsm (rsb_real_t * sfp, rsb_int_t * tnp, rsb_int_t maxr,rsb_time_t maxt, rsb_trans_t transA, const NT * alphap, rsb_coo_idx_tnrhs, rsb_flags_t order, const NT * Bp, rsb_nnz_idx_t ldB, const NT *betap, NT * Cp, rsb_nnz_idx_t ldC) [inline]
template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::tune_spsm_threads (rsb_real_t * sfp = RSBP_NULL, rsb_int_t * tnp =RSBP_NULL, rsb_int_t maxr = 0, rsb_time_t maxt = 0, rsb_trans_t transA= RSB_TRANSPOSITION_N, const NT * alphap = RSBP_NULL, rsb_coo_idx_tnrhs = 1, rsb_flags_t order = RSB_FLAG_WANT_COLUMN_MAJOR_ORDER, constNT * Bp = RSBP_NULL, rsb_nnz_idx_t ldB = 0, const NT * betap =RSBP_NULL, NT * Cp = RSBP_NULL, rsb_nnz_idx_t ldC = 0) const [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::upd_vals (enumrsb_elopf_t elop_flags, const NT & omega) [inline]
template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::upd_vals (enumrsb_elopf_t elop_flags, const NT * omegap) [inline]
Author

NAME

RsbMatrix - Represent a sparse matrix in RSB format by means of librsb .

SYNOPSIS

#include <rsb.hpp>

Public Types

enum RsbSym { IsGen = RSB_FLAG_NOFLAGS, IsHer = RSB_FLAG_HERMITIAN, IsSym = RSB_FLAG_SYMMETRIC, IsTri = RSB_FLAG_TRIANGULAR }

Public Member Functions

RsbMatrix ( rsb_coo_idx_t nrA , rsb_coo_idx_t ncA , const RsbSym sym = IsGen )
RsbMatrix
( rsb_coo_idx_t nrA , const rsb_coo_idx_t * RP , const rsb_coo_idx_t * JA , const NT * VA , const RsbSym sym = IsGen )
RsbMatrix
( const rsb_coo_idx_t * IA , const rsb_coo_idx_t * JA , const NT * VA , rsb_nnz_idx_t nnzA , const rsb_flags_t flagsA = RSB_FLAG_NOFLAGS )
RsbMatrix
( const rsb_char_t * filename , const RsbSym sym = IsGen )
RsbMatrix
( const RsbMatrix & A_Rsb , bool do_trans = false , rsb_flags_t flagsA = RSB_FLAG_NOFLAGS )
RsbMatrix
( RsbMatrix && other )
˜RsbMatrix
( void )
RSBP_RVT RSBP_DEPRECATED Err_t _add
( rsb_coo_idx_t i , rsb_coo_idx_t j , NT val )
RSBP_RVT Err_t close
( void )
RSBP_RVT RSBP_DEPRECATED Err_t _close
( void )
RSBP_RVT RSBP_DEPRECATED Err_t spmv
( rsb_trans_t transA , const NT * alphap , const NT * Xp , rsb_coo_idx_t incX , const NT * betap , NT * Yp , rsb_coo_idx_t incY ) const
RSBP_RVT Err_t spmv
( rsb_trans_t transA , const NT alpha , const NT * Xp , rsb_coo_idx_t incX , const NT beta , NT * Yp , rsb_coo_idx_t incY ) const
RSBP_RVT Err_t spmv
( rsb_trans_t transA , const NT alpha , const NT * Xp , const NT beta , NT * Yp ) const
RSBP_RVT Err_t spmv
(NT * y , const NT * x , bool do_trans = false ) const
RSBP_RVT RSBP_DEPRECATED Err_t spmm
( rsb_trans_t transA , const NT * alphap , rsb_coo_idx_t nrhs , rsb_flags_t order , const NT * Bp , rsb_nnz_idx_t ldB , const NT * betap , NT * Cp , rsb_nnz_idx_t ldC ) const
RSBP_RVT Err_t spmm
( rsb_trans_t transA , const NT alpha , rsb_coo_idx_t nrhs , rsb_flags_t order , const NT * Bp , rsb_nnz_idx_t ldB , const NT beta , NT * Cp , rsb_nnz_idx_t ldC ) const
RSBP_RVT Err_t spmm
( rsb_trans_t transA , const NT alpha , rsb_coo_idx_t nrhs , rsb_flags_t order , const NT * Bp , const NT beta , NT * Cp ) const
RSBP_RVT RSBP_DEPRECATED Err_t spsm
( rsb_trans_t transT , const NT * alphap , rsb_coo_idx_t nrhs , rsb_flags_t order , const NT * betap , const NT * Bp , rsb_nnz_idx_t ldB , NT * Cp , rsb_nnz_idx_t ldC ) const
RSBP_RVT Err_t spsm
( rsb_trans_t transT , const NT alpha , rsb_coo_idx_t nrhs , rsb_flags_t order , const NT beta , const NT * Bp , rsb_nnz_idx_t ldB , NT * Cp , rsb_nnz_idx_t ldC ) const
RSBP_RVT Err_t spsm
( rsb_trans_t transT , const NT alpha , rsb_coo_idx_t nrhs , const NT * Bp , NT * Cp ) const
RSBP_RVT Err_t spsm
(NT * y , const NT * x , rsb_coo_idx_t nrhs , bool do_trans = false ) const
RSBP_RVT Err_t spsm
(NT * y , rsb_coo_idx_t nrhs , bool do_trans = false ) const
RSBP_RVT RSBP_DEPRECATED Err_t spsv
( rsb_trans_t transT , const NT * alphap , const NT * Xp , rsb_coo_idx_t incX , NT * Yp , rsb_coo_idx_t incY ) const
RSBP_RVT Err_t spsv
( rsb_trans_t transT , const NT alpha , const NT * Xp , NT * Yp ) const
RSBP_RVT Err_t spsv
(NT * y , const NT * x , bool do_trans = false ) const
RSBP_RVT Err_t spsv
(NT * y , bool do_trans = false ) const
size_t get_info_size_t
( enum rsb_mif_t mif ) const
rsb_flags_t get_info_rsb_flags_t
( enum rsb_mif_t mif ) const
rsb_blk_idx_t get_info_blk_t
( enum rsb_mif_t mif ) const
rsb_nnz_idx_t get_info_nnz_t
( enum rsb_mif_t mif ) const
rsb_flags_t get_flags_t
( enum rsb_mif_t mif ) const
rsb_type_t get_type_t
( enum rsb_mif_t mif ) const
rsb_coo_idx_t get_info_coo_t
( enum rsb_mif_t mif ) const
size_t _get_index_storage_bytes
( void ) const
size_t _get_storage_bytes
( void ) const
rsb_nnz_idx_t nnz
( void ) const
rsb_blk_idx_t blocks
( void ) const
rsb_coo_idx_t rows
( void ) const
rsb_coo_idx_t cols
( void ) const
RSBP_RVT Err_t get_vals
(NT * VA , const rsb_coo_idx_t * IA , const rsb_coo_idx_t * JA , rsb_nnz_idx_t nnz , rsb_flags_t flags ) const
NT get_val ( const rsb_coo_idx_t i , const rsb_coo_idx_t j , rsb_flags_t flags = RSB_FLAG_NOFLAGS ) const
RSBP_RVT Err_t set_val
( const NT val , const rsb_coo_idx_t i , const rsb_coo_idx_t j , rsb_flags_t flags = RSB_FLAG_NOFLAGS )
RSBP_RVT Err_t set_vals
( const NT * VA , const rsb_coo_idx_t * IA , const rsb_coo_idx_t * JA , rsb_nnz_idx_t nnz , rsb_flags_t flags )
RSBP_RVT Err_t get_vec
(NT * Dp , enum rsb_extff_t flags ) const
RSBP_RVT RSBP_DEPRECATED Err_t get_coo
( rsb_trans_t transA , NT * VA , rsb_coo_idx_t * IA , rsb_coo_idx_t * JA , rsb_flags_t flags ) const
RSBP_RVT RSBP_DEPRECATED Err_t get_csr
( rsb_trans_t transA , NT * VA , rsb_coo_idx_t * RP , rsb_coo_idx_t * JA , rsb_flags_t flags ) const
RSBP_RVT RSBP_DEPRECATED Err_t get_rows_sparse
( rsb_trans_t transA , const NT * alphap , NT * VA , rsb_coo_idx_t * IA , rsb_coo_idx_t * JA , rsb_coo_idx_t frA , rsb_coo_idx_t lrA , rsb_nnz_idx_t * rnzp , rsb_flags_t flags ) const
RSBP_RVT Err_t upd_vals
( enum rsb_elopf_t elop_flags , const NT & omega )
RSBP_RVT Err_t upd_vals
( enum rsb_elopf_t elop_flags , const NT * omegap )
RSBP_RVT Err_t get_nrm
(NT * Np , enum rsb_extff_t flags ) const
rsb_type_t rsbtype
( void ) const
rsb_flags_t rsbflags
( void ) const
rsb_string_t get_info_str
( const char * key ) const
RSBP_RVT Err_t get_info
( enum rsb_mif_t miflags , void * minfop ) const
rsb_string_t _info
( void ) const
RSBP_RVT RSBP_DEPRECATED Err_t tune_spsm_threads
( rsb_real_t * sfp = RSBP_NULL , rsb_int_t * tnp = RSBP_NULL , rsb_int_t maxr =0, rsb_time_t maxt =0, rsb_trans_t transA = RSB_TRANSPOSITION_N , const NT * alphap = RSBP_NULL , rsb_coo_idx_t nrhs =1, rsb_flags_t order = RSB_FLAG_WANT_COLUMN_MAJOR_ORDER , const NT * Bp = RSBP_NULL , rsb_nnz_idx_t ldB =0, const NT * betap = RSBP_NULL , NT * Cp = RSBP_NULL , rsb_nnz_idx_t ldC =0) const
RSBP_RVT RSBP_DEPRECATED Err_t tune_spmm_threads
( rsb_real_t * sfp = RSBP_NULL , rsb_int_t * tnp = RSBP_NULL , rsb_int_t maxr =0, rsb_time_t maxt =0, rsb_trans_t transA = RSB_TRANSPOSITION_N , const NT * alphap = RSBP_NULL , rsb_coo_idx_t nrhs =1, rsb_flags_t order = RSB_FLAG_WANT_COLUMN_MAJOR_ORDER , const NT * Bp = RSBP_NULL , rsb_nnz_idx_t ldB =0, const NT * betap = RSBP_NULL , NT * Cp = RSBP_NULL , rsb_nnz_idx_t ldC =0) const
RSBP_RVT RSBP_DEPRECATED Err_t tune_spmm
( rsb_real_t * sfp , rsb_int_t * tnp , rsb_int_t maxr , rsb_time_t maxt , rsb_trans_t transA , const NT * alphap , rsb_coo_idx_t nrhs , rsb_flags_t order , const NT * Bp , rsb_nnz_idx_t ldB , const NT * betap , NT * Cp , rsb_nnz_idx_t ldC )
RSBP_RVT RSBP_DEPRECATED Err_t tune_spmm
( rsb_real_t * sfp , rsb_int_t * tnp , rsb_int_t maxr , rsb_time_t maxt , rsb_trans_t transA , const NT alpha , rsb_coo_idx_t nrhs , rsb_flags_t order , const NT * Bp , rsb_nnz_idx_t ldB , const NT beta , NT * Cp , rsb_nnz_idx_t ldC )
RSBP_RVT Err_t tune_spmm
( rsb_real_t & sf , rsb_trans_t transA , const NT alpha , rsb_coo_idx_t nrhs , rsb_flags_t order , const NT * Bp , const NT beta , NT * Cp )
RSBP_RVT RSBP_DEPRECATED Err_t tune_spsm
( rsb_real_t * sfp , rsb_int_t * tnp , rsb_int_t maxr , rsb_time_t maxt , rsb_trans_t transA , const NT * alphap , rsb_coo_idx_t nrhs , rsb_flags_t order , const NT * Bp , rsb_nnz_idx_t ldB , const NT * betap , NT * Cp , rsb_nnz_idx_t ldC )
RSBP_RVT Err_t file_save
( const rsb_char_t * filename = RSBP_NULL ) const
RsbMatrix
& operator= ( const RsbMatrix & A_Rsb )
bool _is_complex
( void ) const
bool operator==
( const RsbMatrix & B_Rsb ) const
bool operator!=
( const RsbMatrix & B_Rsb ) const
NT normOne ( void ) const
NT normInf ( void ) const
RSBP_RVT Err_t rndr
( const rsb_char_t * filename = RSBP_NULL , rsb_coo_idx_t pmWidth =512, rsb_coo_idx_t pmHeight =512, rsb_marf_t rflags = RSB_MARF_EPS ) const

Detailed Description

template<RSBP_Scalar_t NT>

class RsbMatrix< NT >"Represent a sparse matrix in RSB format by means of librsb .

Manage construction, destruction, and std::move of numerical matrices.
Most of the member functions here translate directly to a single function call to librsb ( rsb.h), and pass the parameters as they are, so the error checking is done by librsb .
While most of
librsb C functions use void* pointers instead of numerical data, RsbMatrix is templated by a type parameter. This introduces type safety at compile time.

Users of member functions can choose among several overloads. So in additional to the more direct overloads passing e.g. $ alpha $ and $ beta $ by reference, here a user can pass them by value.

Parameters

NT the numerical type, at least for the four canonical ones ( float , double , std::complex<float> , std::complex<double> ); see matrix_supported_numerical_types_section and rsb_type_t for more.

Note

Default error propagation is by exception throw for all constructors and most member functions.
Functions declared to return Err_t can be specialized in rsb_err_t so not to throw exceptions, but to return an error code instead.
Exceptions thrown by member functions (not constructors) can be deactivated at build time by defining RSBP_NOTHROW before including < rsb.hpp > .

One may turn on return error value as default at build time by defining RSBP_WANT_REV=1 .

Warning

The error model is work in progress and subject to change.

Todo

While the rsb.h interface is stable, the rsb.hpp interface is neither stable, nor complete: early users’ feedback is very welcome.
Shall all $ alpha $ and $ beta $ be passed by values only? This is natural and fits C++. But rsb_tune_spmm() / rsb_tune_spsm() have more parameters with a nullptr-like ‘default’ value: how to deal with them in RsbMatrix::tune_spmm() / RsbMatrix::tune_spsm() without introducing too many special cases ?
Similarly for the consistency of RsbMatrix::get_flags_t(), RsbMatrix::get_type_t(), RsbMatrix::get_info_coo_t(), and similar: shall they throw an exception if given a flag not matching the return type ?
Or maybe keeping only RsbMatrix::get_info(), with per-type reference overloads ?
Furthermore: if working under C++20, shall RsbLib avoid pointer-based interfaces completely (using
std::span only)?

Member Enumeration Documentation

template<RSBP_Scalar_t NT> enum RsbMatrix::RsbSym

Matrix structure: either general, symmetric, hermitian, or triangular (also see rsb_flags_t ).

Enumerator

IsGen

General matrix, no triangle structure or symmetry assumed.

IsHer

Hermitian ( $ A == AˆH $). Please pass only lower/upper triangle.

IsSym

Symmetric ( $ A == AˆT $). Please pass only lower/upper triangle.

IsTri

Triangular (required for spsv/ spsm).

Constructor & Destructor Documentation

template<RSBP_Scalar_t NT> RsbMatrix< NT >::RsbMatrix (rsb_coo_idx_t nrA,rsb_coo_idx_t ncA, const RsbSym sym = IsGen) [inline]

Begin assembling a sparse matrix of given dimensions and type.
Then you can use set_val() or set_vals() repeatedly to populate the matrix.
After populating the matrix, use close() to terminate its assembly.
Example snip from examples/assemble.cpp :

const rsb_coo_idx_t nrA { 4 }, ncA { 4 };
RsbMatrix<double> mtx(nrA,ncA); // begin matrix assembly

// insert elements of a tridiagonal matrix, one by one
for (auto i = 0; i < nrA; ++i )
for (auto j = i-1; j <= i+1; ++j )
if ( i >= 0 && i < nrA )
if ( j >= 0 && j < ncA )
mtx.set_val((i+1)*100+(j+1),i,j); // add entry

mtx.close(); // finish matrix assembly
assert(mtx.nnz() == 3 * nrA - 2);

Note

Directly based on C function rsb_mtx_alloc_from_coo_begin() from < rsb.h > .

template<RSBP_Scalar_t NT> RsbMatrix< NT >::RsbMatrix (rsb_coo_idx_t nrA,const rsb_coo_idx_t * RP, const rsb_coo_idx_t * JA, const NT * VA,const RsbSym sym = IsGen) [inline]

Assemble a sparse matrix given CSR input.
Note

Directly based on C function rsb_mtx_alloc_from_csr_const() from < rsb.h > .

template<RSBP_Scalar_t NT> RsbMatrix< NT >::RsbMatrix (const rsb_coo_idx_t* IA, const rsb_coo_idx_t * JA, const NT * VA, rsb_nnz_idx_t nnzA,const rsb_flags_t flagsA = RSB_FLAG_NOFLAGS) [inline]

Assemble a sparse matrix given COO input.
Example snip from examples/misc.cpp :

const rsb_nnz_idx_t nnzA { 7 };
const rsb_coo_idx_t nrA { 6 }, ncA { 6 }, nrhs { 1 };
const std::vector<rsb_coo_idx_t> IA {0,1,2,3,4,5,1}, JA {0,1,2,3,4,5,0};
const std::vector<double> VA {1,1,1,1,1,1,2}, X(ncA,1);
RsbMatrix<double> mtx(IA.data(),JA.data(),VA.data(),nnzA);

Note

Directly based on C function rsb_mtx_alloc_from_coo_const() from < rsb.h > .

template<RSBP_Scalar_t NT> RsbMatrix< NT >::RsbMatrix (const rsb_char_t *filename, const RsbSym sym = IsGen) [inline]

Assemble a sparse matrix given filename input.
Note

Directly based on C function rsb_file_mtx_load() from < rsb.h > .

template<RSBP_Scalar_t NT> RsbMatrix< NT >::RsbMatrix (const RsbMatrix< NT> & A_Rsb, bool do_trans = false, rsb_flags_t flagsA =RSB_FLAG_NOFLAGS) [inline]

Copy a sparse matrix given example input.
Can either clone it, or transpose it or change flags (structure) in the process.

template<RSBP_Scalar_t NT> RsbMatrix< NT >::RsbMatrix (RsbMatrix< NT > &&other) [inline]

Move constructor.
The moved matrix object will be invalid afterwards.
Example snip from examples/misc.cpp :

assert( mtx1.nnz() == nnzA );
RsbMatrix<double> mtx3 { std::move(mtx1) };
assert( mtx3.nnz() == nnzA );

template<RSBP_Scalar_t NT> RsbMatrix< NT >::˜RsbMatrix (void) [inline]

Destructor.
Frees matrix object memory.
Note

Directly based on C function rsb_mtx_free() from < rsb.h > .

Member Function Documentation

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::_add (rsb_coo_idx_t i, rsb_coo_idx_t j, NT val) [inline]

Deprecated

Use set_val() and set_vals() instead.

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::_close (void) [inline]

Deprecated

Use close() instead.

template<RSBP_Scalar_t NT> size_t RsbMatrix< NT >::_get_index_storage_bytes(void) const [inline]

Warning

The name of this member function is expected to change.

template<RSBP_Scalar_t NT> size_t RsbMatrix< NT >::_get_storage_bytes(void) const [inline]

Warning

The name of this member function is expected to change.

template<RSBP_Scalar_t NT> rsb_string_t RsbMatrix< NT >::_info (void) const[inline]

Warning

The name of this member function is expected to change.

template<RSBP_Scalar_t NT> bool RsbMatrix< NT >::_is_complex (void) const[inline]

Warning

The name of this member function is expected to change.

template<RSBP_Scalar_t NT> rsb_blk_idx_t RsbMatrix< NT >::blocks (void)const [inline]

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::close (void)[inline]

Terminate assembly of a previously started and populated matrix.
Shall be called once.

Example snip from examples/assemble.cpp :

const rsb_coo_idx_t nrA { 4 }, ncA { 4 };
RsbMatrix<double> mtx(nrA,ncA); // begin matrix assembly

// insert elements of a tridiagonal matrix, one by one
for (auto i = 0; i < nrA; ++i )
for (auto j = i-1; j <= i+1; ++j )
if ( i >= 0 && i < nrA )
if ( j >= 0 && j < ncA )
mtx.set_val((i+1)*100+(j+1),i,j); // add entry

mtx.close(); // finish matrix assembly
assert(mtx.nnz() == 3 * nrA - 2);

See also

RsbMatrix::RsbMatrix ( rsb_coo_idx_t nrA , rsb_coo_idx_t ncA , const RsbSym sym = IsGen );

Note

Directly based on C function rsb_mtx_alloc_from_coo_end() from < rsb.h > .

template<RSBP_Scalar_t NT> rsb_coo_idx_t RsbMatrix< NT >::cols (void) const[inline]

Example snip from examples/twonnz.cpp :

std::cout << "# Matrix sized " << mtx.rows() << "x" << mtx.cols() << ", " << nnzA << " nnz built in " << dt << " s and occupies " << mtxocc << " bytes " << std::endl;

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::file_save (constrsb_char_t * filename = RSBP_NULL) const [inline]

Note

Directly based on C function rsb_file_mtx_save() from < rsb.h > .

Example snip from examples/misc.cpp :

mtx.file_save(); // print to stdout

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::get_coo (rsb_trans_t transA, NT * VA, rsb_coo_idx_t * IA,rsb_coo_idx_t * JA, rsb_flags_t flags) const [inline]

Note

Directly based on C function rsb_mtx_get_coo() from < rsb.h > .

Warning

Only transA=RSB_TRANSPOSITION_N currently supported.

Deprecated

This function overload is deprecated: use the corresponding overload without the C-style pointers-to-scalar, or request one if non existing yet.

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::get_csr (rsb_trans_t transA, NT * VA, rsb_coo_idx_t * RP,rsb_coo_idx_t * JA, rsb_flags_t flags) const [inline]

Note

Directly based on C function rsb_mtx_get_csr() from < rsb.h > .

Warning

Only transA=RSB_TRANSPOSITION_N currently supported.

Deprecated

This function overload is deprecated: use the corresponding overload without the C-style pointers-to-scalar, or request one if non existing yet.

template<RSBP_Scalar_t NT> rsb_flags_t RsbMatrix< NT >::get_flags_t (enumrsb_mif_t mif) const [inline]

Note

Directly based on C function rsb_mtx_get_info() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::get_info (enumrsb_mif_t miflags, void * minfop) const [inline]

Note

Directly based on C function rsb_mtx_get_info() from < rsb.h > .

template<RSBP_Scalar_t NT> rsb_blk_idx_t RsbMatrix< NT >::get_info_blk_t(enum rsb_mif_t mif) const [inline]

Note

Directly based on C function rsb_mtx_get_info() from < rsb.h > .

template<RSBP_Scalar_t NT> rsb_coo_idx_t RsbMatrix< NT >::get_info_coo_t(enum rsb_mif_t mif) const [inline]

Note

Directly based on C function rsb_mtx_get_info() from < rsb.h > .

template<RSBP_Scalar_t NT> rsb_nnz_idx_t RsbMatrix< NT >::get_info_nnz_t(enum rsb_mif_t mif) const [inline]

Note

Directly based on C function rsb_mtx_get_info() from < rsb.h > .

template<RSBP_Scalar_t NT> rsb_flags_t RsbMatrix< NT>::get_info_rsb_flags_t (enum rsb_mif_t mif) const [inline]

Note

Directly based on C function rsb_mtx_get_info() from < rsb.h > .

template<RSBP_Scalar_t NT> size_t RsbMatrix< NT >::get_info_size_t (enumrsb_mif_t mif) const [inline]

Note

Directly based on C function rsb_mtx_get_info() from < rsb.h > .

template<RSBP_Scalar_t NT> rsb_string_t RsbMatrix< NT >::get_info_str(const char * key) const [inline]

Note

Directly based on C function rsb_mtx_get_info_str() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::get_nrm (NT *Np, enum rsb_extff_t flags) const [inline]

Note

Directly based on C function rsb_mtx_get_nrm() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::get_rows_sparse (rsb_trans_t transA, const NT * alphap, NT * VA,rsb_coo_idx_t * IA, rsb_coo_idx_t * JA, rsb_coo_idx_t frA,rsb_coo_idx_t lrA, rsb_nnz_idx_t * rnzp, rsb_flags_t flags) const[inline]

Note

Directly based on C function rsb_mtx_get_rows_sparse() from < rsb.h > .

Deprecated

This function overload is deprecated: use the corresponding overload without the C-style pointers-to-scalar, or request one if non existing yet.

template<RSBP_Scalar_t NT> rsb_type_t RsbMatrix< NT >::get_type_t (enumrsb_mif_t mif) const [inline]

Note

Directly based on C function rsb_mtx_get_info() from < rsb.h > .

template<RSBP_Scalar_t NT> NT RsbMatrix< NT >::get_val (const rsb_coo_idx_ti, const rsb_coo_idx_t j, rsb_flags_t flags = RSB_FLAG_NOFLAGS) const[inline]

Note

Directly based on C function rsb_mtx_get_vals() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::get_vals (NT *VA, const rsb_coo_idx_t * IA, const rsb_coo_idx_t * JA, rsb_nnz_idx_tnnz, rsb_flags_t flags) const [inline]

Note

Directly based on C function rsb_mtx_get_vals() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::get_vec (NT *Dp, enum rsb_extff_t flags) const [inline]

Note

Directly based on C function rsb_mtx_get_vec() from < rsb.h > .

template<RSBP_Scalar_t NT> rsb_nnz_idx_t RsbMatrix< NT >::nnz (void) const[inline]

template<RSBP_Scalar_t NT> NT RsbMatrix< NT >::normInf (void) const[inline]

template<RSBP_Scalar_t NT> NT RsbMatrix< NT >::normOne (void) const[inline]

template<RSBP_Scalar_t NT> bool RsbMatrix< NT >::operator!= (constRsbMatrix< NT > & B_Rsb) const [inline]

Example snip from examples/misc.cpp :

assert( mtx1 == mtx2 );
assert( !(mtx1 != mtx2) );

See also

RsbMatrix::operator==(const RsbMatrix & B_Rsb) const ;

template<RSBP_Scalar_t NT> RsbMatrix & RsbMatrix< NT >::operator= (constRsbMatrix< NT > & A_Rsb) [inline]

A copy constructor. Will clone the input matrix contents.

template<RSBP_Scalar_t NT> bool RsbMatrix< NT >::operator== (constRsbMatrix< NT > & B_Rsb) const [inline]

Deep comparison: compare if the two matrices have same dimensions, nonzeroes count, nonzeroes pattern and value. Meant for very sporadic use. Inefficient: it can involve matrices copying.

Example snip from examples/misc.cpp :

assert( mtx1 == mtx2 );
assert( !(mtx1 != mtx2) );

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::rndr (constrsb_char_t * filename = RSBP_NULL, rsb_coo_idx_t pmWidth = 512,rsb_coo_idx_t pmHeight = 512, rsb_marf_t rflags = RSB_MARF_EPS) const[inline]

Note

Directly based on C function rsb_mtx_rndr() from < rsb.h > .

Example snip from examples/render.cpp :

mtx.rndr(psfilename.c_str());

template<RSBP_Scalar_t NT> rsb_coo_idx_t RsbMatrix< NT >::rows (void) const[inline]

Example snip from examples/twonnz.cpp :

std::cout << "# Matrix sized " << mtx.rows() << "x" << mtx.cols() << ", " << nnzA << " nnz built in " << dt << " s and occupies " << mtxocc << " bytes " << std::endl;

template<RSBP_Scalar_t NT> rsb_flags_t RsbMatrix< NT >::rsbflags (void)const [inline]

template<RSBP_Scalar_t NT> rsb_type_t RsbMatrix< NT >::rsbtype (void) const[inline]

Example snip from examples/twonnz.cpp :

std::cout << "# type=" << mtx.rsbtype() << " nt=1," << rnt << " n=" << n << " nrhs=" << nrhs << " order=" << oc << " alpha=" << alpha << " beta=" << beta << " dt=" << dta[0] << ".." << dta[1] << " spmm-scalability=" << dta[0]/dta[1] << " nnz/s=" << nnzA/dta[0] << ".." << nnzA/dta[1] << " flops=" << flops_c/dta[0] << ".." << flops_c/dta[1] << " occ.=" << opocc << " " << std::endl;

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::set_val (constNT val, const rsb_coo_idx_t i, const rsb_coo_idx_t j, rsb_flags_t flags= RSB_FLAG_NOFLAGS) [inline]

Note

Directly based on C function rsb_mtx_set_vals() from < rsb.h > .

Example snip from examples/assemble.cpp :

const rsb_coo_idx_t nrA { 4 }, ncA { 4 };
RsbMatrix<double> mtx(nrA,ncA); // begin matrix assembly

// insert elements of a tridiagonal matrix, one by one
for (auto i = 0; i < nrA; ++i )
for (auto j = i-1; j <= i+1; ++j )
if ( i >= 0 && i < nrA )
if ( j >= 0 && j < ncA )
mtx.set_val((i+1)*100+(j+1),i,j); // add entry

mtx.close(); // finish matrix assembly
assert(mtx.nnz() == 3 * nrA - 2);

See also

set_vals() .

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::set_vals (constNT * VA, const rsb_coo_idx_t * IA, const rsb_coo_idx_t * JA,rsb_nnz_idx_t nnz, rsb_flags_t flags) [inline]

Add a single entry during the assembly of a matrix created empty.
Use close() to terminate matrix assembly.

Example snip from examples/assemble.cpp :

const rsb_coo_idx_t nrA { 4 }, ncA { 4 };
RsbMatrix<double> mtx(nrA,ncA); // begin matrix assembly

// insert elements of a tridiagonal matrix, one by one
for (auto i = 0; i < nrA; ++i )
for (auto j = i-1; j <= i+1; ++j )
if ( i >= 0 && i < nrA )
if ( j >= 0 && j < ncA )
mtx.set_val((i+1)*100+(j+1),i,j); // add entry

mtx.close(); // finish matrix assembly
assert(mtx.nnz() == 3 * nrA - 2);

Note

Directly based on C function rsb_mtx_set_vals() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::spmm (rsb_trans_t transA, const NT * alphap, rsb_coo_idx_t nrhs,rsb_flags_t order, const NT * Bp, rsb_nnz_idx_t ldB, const NT * betap,NT * Cp, rsb_nnz_idx_t ldC) const [inline]

Note

Directly based on C function rsb_spmm() from < rsb.h > .

Deprecated

This function overload is deprecated: use the corresponding overload without the C-style pointers-to-scalar, or request one if non existing yet.

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spmm(rsb_trans_t transA, const NT alpha, rsb_coo_idx_t nrhs, rsb_flags_torder, const NT * Bp, const NT beta, NT * Cp) const [inline]

Note

Indirectly based on C function rsb_spmm() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spmm(rsb_trans_t transA, const NT alpha, rsb_coo_idx_t nrhs, rsb_flags_torder, const NT * Bp, rsb_nnz_idx_t ldB, const NT beta, NT * Cp,rsb_nnz_idx_t ldC) const [inline]

Example snip from examples/bench.cpp :

mtx.spmm(transA,&alpha,nrhs,order,B.data(),ldB,&beta,C.data(),ldC);

Note

Indirectly based on C function rsb_spmm() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spmv (NT * y,const NT * x, bool do_trans = false) const [inline]

Note

Indirectly based on C function rsb_spmv() from < rsb.h > .

Example snip from examples/span.cpp :

mtx.tune_spmm(nullptr,&tn,0,0.0,RSB_TRANSPOSITION_N,alpha,nrhs,RSB_FLAG_WANT_COLUMN_MAJOR_ORDER,X,ncA,beta,Y,nrA);
mtx.spmv(RSB_TRANSPOSITION_N, alpha, X, beta, Y);

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::spmv (rsb_trans_t transA, const NT * alphap, const NT * Xp,rsb_coo_idx_t incX, const NT * betap, NT * Yp, rsb_coo_idx_t incY)const [inline]

Note

Directly based on C function rsb_spmv() from < rsb.h > .

Deprecated

This function overload is deprecated: use the corresponding overload without the C-style pointers-to-scalar, or request one if non existing yet.

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spmv(rsb_trans_t transA, const NT alpha, const NT * Xp, const NT beta, NT *Yp) const [inline]

Note

Indirectly based on C function rsb_spmv() from < rsb.h > .

Example snip from examples/span.cpp :

mtx.tune_spmm(nullptr,&tn,0,0.0,RSB_TRANSPOSITION_N,alpha,nrhs,RSB_FLAG_WANT_COLUMN_MAJOR_ORDER,X,ncA,beta,Y,nrA);
mtx.spmv(RSB_TRANSPOSITION_N, alpha, X, beta, Y);

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spmv(rsb_trans_t transA, const NT alpha, const NT * Xp, rsb_coo_idx_t incX,const NT beta, NT * Yp, rsb_coo_idx_t incY) const [inline]

Note

Indirectly based on C function rsb_spmv() from < rsb.h > .

Example snip from examples/span.cpp :

mtx.tune_spmm(nullptr,&tn,0,0.0,RSB_TRANSPOSITION_N,alpha,nrhs,RSB_FLAG_WANT_COLUMN_MAJOR_ORDER,X,ncA,beta,Y,nrA);
mtx.spmv(RSB_TRANSPOSITION_N, alpha, X, beta, Y);

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsm (NT * y,const NT * x, rsb_coo_idx_t nrhs, bool do_trans = false) const [inline]

Note

Indirectly based on C function rsb_spsm() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsm (NT * y,rsb_coo_idx_t nrhs, bool do_trans = false) const [inline]

Note

Indirectly based on C function rsb_spsm() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::spsm (rsb_trans_t transT, const NT * alphap, rsb_coo_idx_t nrhs,rsb_flags_t order, const NT * betap, const NT * Bp, rsb_nnz_idx_t ldB,NT * Cp, rsb_nnz_idx_t ldC) const [inline]

Note

Directly based on C function rsb_spsm() from < rsb.h > .

Deprecated

This function overload is deprecated: use the corresponding overload without the C-style pointers-to-scalar, or request one if non existing yet.

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsm(rsb_trans_t transT, const NT alpha, rsb_coo_idx_t nrhs, const NT * Bp,NT * Cp) const [inline]

Note

Indirectly based on C function rsb_spsm() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsm(rsb_trans_t transT, const NT alpha, rsb_coo_idx_t nrhs, rsb_flags_torder, const NT beta, const NT * Bp, rsb_nnz_idx_t ldB, NT * Cp,rsb_nnz_idx_t ldC) const [inline]

Note

Indirectly based on C function rsb_spsm() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsv (NT * y,bool do_trans = false) const [inline]

Note

Indirectly based on C function rsb_spsv() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsv (NT * y,const NT * x, bool do_trans = false) const [inline]

Note

Indirectly based on C function rsb_spsv() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::spsv (rsb_trans_t transT, const NT * alphap, const NT * Xp,rsb_coo_idx_t incX, NT * Yp, rsb_coo_idx_t incY) const [inline]

Note

Directly based on C function rsb_spsv() from < rsb.h > .

Deprecated

This function overload is deprecated: use the corresponding overload without the C-style pointers-to-scalar, or request one if non existing yet.

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::spsv(rsb_trans_t transT, const NT alpha, const NT * Xp, NT * Yp) const[inline]

Note

Indirectly based on C function rsb_spsv() from < rsb.h > .

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::tune_spmm(rsb_real_t & sf, rsb_trans_t transA, const NT alpha, rsb_coo_idx_tnrhs, rsb_flags_t order, const NT * Bp, const NT beta, NT * Cp)[inline]

Note

Indirectly based on C function rsb_tune_spmm() from < rsb.h > .

Example snip from examples/autotune.cpp :

tt = -rsb_time();
mtx.tune_spmm(&sf,&tn,maxr,tmax,transA,&alpha,nrhs,order,B.data(),ldB,&beta,C.data(),ldC);
tt += rsb_time();

auto nnsmA {mtx.blocks()};
std::cout << "Tuning took " << tt << " s ( " << tt / dt << " ops ) and changed " << nsmA << " to " << nnsmA << " blocks" << std::endl;

mtx.spmm(transA,&alpha,nrhs,order,B.data(),ldB,&beta,C.data(),ldC); // caches warmup

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::tune_spmm (rsb_real_t * sfp, rsb_int_t * tnp, rsb_int_t maxr,rsb_time_t maxt, rsb_trans_t transA, const NT * alphap, rsb_coo_idx_tnrhs, rsb_flags_t order, const NT * Bp, rsb_nnz_idx_t ldB, const NT *betap, NT * Cp, rsb_nnz_idx_t ldC) [inline]

Note

Directly based on C function rsb_tune_spmm() from < rsb.h > .

Example snip from examples/autotune.cpp :

tt = -rsb_time();
mtx.tune_spmm(&sf,&tn,maxr,tmax,transA,&alpha,nrhs,order,B.data(),ldB,&beta,C.data(),ldC);
tt += rsb_time();

auto nnsmA {mtx.blocks()};
std::cout << "Tuning took " << tt << " s ( " << tt / dt << " ops ) and changed " << nsmA << " to " << nnsmA << " blocks" << std::endl;

mtx.spmm(transA,&alpha,nrhs,order,B.data(),ldB,&beta,C.data(),ldC); // caches warmup

Deprecated

This function overload is deprecated: use the corresponding overload without the C-style pointers-to-scalar, or request one if non existing yet.

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::tune_spmm (rsb_real_t * sfp, rsb_int_t * tnp, rsb_int_t maxr,rsb_time_t maxt, rsb_trans_t transA, const NT alpha, rsb_coo_idx_tnrhs, rsb_flags_t order, const NT * Bp, rsb_nnz_idx_t ldB, const NTbeta, NT * Cp, rsb_nnz_idx_t ldC) [inline]

Note

Indirectly based on C function rsb_tune_spmm() from < rsb.h > .

Example snip from examples/autotune.cpp :

tt = -rsb_time();
mtx.tune_spmm(&sf,&tn,maxr,tmax,transA,&alpha,nrhs,order,B.data(),ldB,&beta,C.data(),ldC);
tt += rsb_time();

auto nnsmA {mtx.blocks()};
std::cout << "Tuning took " << tt << " s ( " << tt / dt << " ops ) and changed " << nsmA << " to " << nnsmA << " blocks" << std::endl;

mtx.spmm(transA,&alpha,nrhs,order,B.data(),ldB,&beta,C.data(),ldC); // caches warmup

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::tune_spmm_threads (rsb_real_t * sfp = RSBP_NULL, rsb_int_t * tnp =RSBP_NULL, rsb_int_t maxr = 0, rsb_time_t maxt = 0, rsb_trans_t transA= RSB_TRANSPOSITION_N, const NT * alphap = RSBP_NULL, rsb_coo_idx_tnrhs = 1, rsb_flags_t order = RSB_FLAG_WANT_COLUMN_MAJOR_ORDER, constNT * Bp = RSBP_NULL, rsb_nnz_idx_t ldB = 0, const NT * betap =RSBP_NULL, NT * Cp = RSBP_NULL, rsb_nnz_idx_t ldC = 0) const [inline]

Note

Directly based on C function rsb_tune_spmm() from < rsb.h > .

Example snip from examples/autotune.cpp :

tt = -rsb_time();
mtx.tune_spmm(&sf,&tn,maxr,tmax,transA,&alpha,nrhs,order,B.data(),ldB,&beta,C.data(),ldC);
tt += rsb_time();

auto nnsmA {mtx.blocks()};
std::cout << "Tuning took " << tt << " s ( " << tt / dt << " ops ) and changed " << nsmA << " to " << nnsmA << " blocks" << std::endl;

mtx.spmm(transA,&alpha,nrhs,order,B.data(),ldB,&beta,C.data(),ldC); // caches warmup

Deprecated

This function overload is deprecated: use the corresponding overload without the C-style pointers-to-scalar, or request one if non existing yet.

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::tune_spsm (rsb_real_t * sfp, rsb_int_t * tnp, rsb_int_t maxr,rsb_time_t maxt, rsb_trans_t transA, const NT * alphap, rsb_coo_idx_tnrhs, rsb_flags_t order, const NT * Bp, rsb_nnz_idx_t ldB, const NT *betap, NT * Cp, rsb_nnz_idx_t ldC) [inline]

Note

Directly based on C function rsb_tune_spsm() from < rsb.h > .

Deprecated

This function overload is deprecated: use the corresponding overload without the C-style pointers-to-scalar, or request one if non existing yet.

template<RSBP_Scalar_t NT> RSBP_RVT RSBP_DEPRECATED Err_t RsbMatrix< NT>::tune_spsm_threads (rsb_real_t * sfp = RSBP_NULL, rsb_int_t * tnp =RSBP_NULL, rsb_int_t maxr = 0, rsb_time_t maxt = 0, rsb_trans_t transA= RSB_TRANSPOSITION_N, const NT * alphap = RSBP_NULL, rsb_coo_idx_tnrhs = 1, rsb_flags_t order = RSB_FLAG_WANT_COLUMN_MAJOR_ORDER, constNT * Bp = RSBP_NULL, rsb_nnz_idx_t ldB = 0, const NT * betap =RSBP_NULL, NT * Cp = RSBP_NULL, rsb_nnz_idx_t ldC = 0) const [inline]

Note

Directly based on C function rsb_tune_spsm() from < rsb.h > .

Deprecated

This function overload is deprecated: use the corresponding overload without the C-style pointers-to-scalar, or request one if non existing yet.

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::upd_vals (enumrsb_elopf_t elop_flags, const NT & omega) [inline]

Note

Directly based on C function rsb_mtx_upd_vals() from < rsb.h > .

Example snip from examples/bench.cpp :

mtx.upd_vals(RSB_ELOPF_POW,nt_t{0.0}); // set matrix values to ones

template<RSBP_Scalar_t NT> RSBP_RVT Err_t RsbMatrix< NT >::upd_vals (enumrsb_elopf_t elop_flags, const NT * omegap) [inline]

Note

Directly based on C function rsb_mtx_upd_vals() from < rsb.h > .

Deprecated

This function overload is deprecated: use the corresponding overload without the C-style pointers-to-scalar, or request one if non existing yet.

Author

Generated automatically by Doxygen for librsb from the source code.