Mentions légales du service

Skip to content
Snippets Groups Projects
Commit 3d3d8d8d authored by hhakim's avatar hhakim
Browse files

Add member functions to MatDense: eq_cols, eq_rows, sum_col, sum_row,...

Add member functions to MatDense: eq_cols, eq_rows, sum_col, sum_row, best_low_rank, row_nonzero_inds, col_nonzero_inds, submatrix.

NOTE: it needs eigen >= 3.4 to compile (slicing support is needed by submatrix()).
parent 8a6d658d
Branches
Tags
No related merge requests found
...@@ -384,6 +384,14 @@ namespace Faust ...@@ -384,6 +384,14 @@ namespace Faust
MatDense<FPP,Cpu>* get_cols(faust_unsigned_int start_col_id, faust_unsigned_int num_cols) const; MatDense<FPP,Cpu>* get_cols(faust_unsigned_int start_col_id, faust_unsigned_int num_cols) const;
MatDense<FPP,Cpu>* get_cols(faust_unsigned_int* col_ids, faust_unsigned_int n) const; MatDense<FPP,Cpu>* get_cols(faust_unsigned_int* col_ids, faust_unsigned_int n) const;
MatDense<FPP,Cpu>* get_cols(std::vector<int> col_ids) const; MatDense<FPP,Cpu>* get_cols(std::vector<int> col_ids) const;
/** \brief Returns true if this[:,id_this] == other[:, id_other] at the specified precision. */
bool eq_cols(const MatDense<FPP, Cpu> & other, faust_unsigned_int id_this, faust_unsigned_int id_other, const Real<FPP>& precision) const;
/** \brief Returns true if this[id_this, :] == other[id_other, :] at the specified precision. */
bool eq_rows(const MatDense<FPP, Cpu> & other, faust_unsigned_int id_this, faust_unsigned_int id_other, const Real<FPP>& precision) const;
/** \brief Returns the sum of the column of index id. */
FPP sum_col(faust_unsigned_int id) const;
/** \brief Returns the sum of the row of index id. */
FPP sum_row(faust_unsigned_int id) const;
void delete_col(int offset); void delete_col(int offset);
void delete_row(int offset); void delete_row(int offset);
...@@ -468,6 +476,22 @@ namespace Faust ...@@ -468,6 +476,22 @@ namespace Faust
*/ */
std::list<std::pair<int,int>> nonzeros_indices() const; std::list<std::pair<int,int>> nonzeros_indices() const;
/**
* \brief Returns the best low rank approximation this = bestX * bestY using the svd.
*/
void best_low_rank(const int &r, MatDense<FPP,Cpu> &bestX, MatDense<FPP, Cpu> &bestY) const;
/**
* Returns the vector of the indices of the nonzeros of the row of index row_id.
*/
std::vector<int> row_nonzero_inds(faust_unsigned_int row_id) const;
/**
* Returns the vector of the indices of the nonzeros of the column of index col_id.
*/
std::vector<int> col_nonzero_inds(faust_unsigned_int col_id) const;
/**
* Returns in submat the submatrix defined by the row indices row_ids and the column indices col_ids.
*/
void submatrix(const std::vector<int> &row_ids, const std::vector<int> &col_ids, MatDense<FPP, Cpu> & submat) const;
private: private:
Eigen::Matrix<FPP, Eigen::Dynamic, Eigen::Dynamic> mat; Eigen::Matrix<FPP, Eigen::Dynamic, Eigen::Dynamic> mat;
bool isZeros; bool isZeros;
......
...@@ -53,7 +53,7 @@ ...@@ -53,7 +53,7 @@
#include "faust_exception.h" #include "faust_exception.h"
#include "faust_constant.h" #include "faust_constant.h"
#include <cassert> #include <cassert>
#include <Eigen/SVD>
namespace Faust namespace Faust
{ {
...@@ -202,7 +202,7 @@ namespace Faust ...@@ -202,7 +202,7 @@ namespace Faust
std::list<std::pair<int,int>> nz_inds; std::list<std::pair<int,int>> nz_inds;
if(this->is_identity) if(this->is_identity)
#ifdef _MSC_VER #ifdef _MSC_VER
for(int i=0;i<min(this->dim1, this->dim2);i++) // VS14 strange issue with std::min //C2589 or C2059 for(int i=0;i<min(this->dim1, this->dim2);i++) // VS14 strange issue with std::min //C2589 or C2059 // TODO/ fix rather with -DNOMINMAX
#else #else
for(int i=0;i<std::min(this->dim1, this->dim2);i++) for(int i=0;i<std::min(this->dim1, this->dim2);i++)
#endif #endif
...@@ -1221,6 +1221,88 @@ MatDense<FPP,Cpu>* MatDense<FPP,Cpu>::get_rows(faust_unsigned_int* row_ids, faus ...@@ -1221,6 +1221,88 @@ MatDense<FPP,Cpu>* MatDense<FPP,Cpu>::get_rows(faust_unsigned_int* row_ids, faus
return rows; return rows;
} }
template<typename FPP>
bool MatDense<FPP,Cpu>::eq_cols(const MatDense<FPP, Cpu> & other, faust_unsigned_int id_this, faust_unsigned_int id_other, const Real<FPP>& precision) const
{
if(this->getNbRow() != other.getNbRow()) return false;
for(int i=0;i<this->getNbRow();i++)
{
if(std::abs((*this)(i, id_this)-other(i, id_other)) > precision)
return false;
}
return true;
}
template<typename FPP>
bool MatDense<FPP,Cpu>::eq_rows(const MatDense<FPP, Cpu> & other, faust_unsigned_int id_this, faust_unsigned_int id_other, const Real<FPP>& precision) const
{
if(this->getNbCol() != other.getNbCol()) return false;
for(int j=0;j<this->getNbCol();j++)
{
if(std::abs((*this)(id_this,j)-other(id_other,j)) > precision)
return false;
}
return true;
}
template<typename FPP>
void MatDense<FPP, Cpu>::best_low_rank(const int &r, MatDense<FPP,Cpu> &bestX, MatDense<FPP, Cpu> &bestY) const
{
Eigen::JacobiSVD<Eigen::Matrix<FPP, Eigen::Dynamic, Eigen::Dynamic>> svd(this->mat, Eigen::ComputeThinU | Eigen::ComputeThinV);
bestX.mat = svd.matrixU().block(0, 0, this->getNbRow(), r) * svd.singularValues().block(0, 0, r, r).asDiagonal();
bestY.mat = svd.matrixV().block(0, 0, this->getNbCol(), r) * svd.singularValues().block(0, 0, r, r).asDiagonal();
}
template<typename FPP>
std::vector<int> MatDense<FPP, Cpu>::col_nonzero_inds(faust_unsigned_int col_id) const
{
std::vector<int> ids;
auto data_ptr = this->getData()+col_id*this->getNbRow();
for(int i=0;i<this->getNbRow();i++)
{
if(data_ptr[i]) ids.push_back(i);
}
return ids;
}
template<typename FPP>
std::vector<int> MatDense<FPP, Cpu>::row_nonzero_inds(faust_unsigned_int row_id) const
{
std::vector<int> ids;
auto data_ptr = this->getData();
for(int j=0;j<this->getNbRow();j++)
if(data_ptr[j*this->getNbRow()+row_id]) ids.push_back(j);
return ids;
}
template<typename FPP>
void MatDense<FPP, Cpu>::submatrix(const std::vector<int> &row_ids, const std::vector<int> &col_ids, MatDense<FPP, Cpu> & submat) const
{
if(this->dim1 != row_ids.size() || this->dim2 != col_ids.size())
submat.resize(row_ids.size(), col_ids.size());
submat.mat = mat(row_ids, col_ids);
std::cout << submat.mat << std::endl;
}
template<typename FPP>
FPP MatDense<FPP,Cpu>::sum_col(faust_unsigned_int id) const
{
FPP sum = (FPP)0;
for(int i=0;i<this->getNbRow();i++)
sum += (*this)(i, id);
return sum;
}
template<typename FPP>
FPP MatDense<FPP,Cpu>::sum_row(faust_unsigned_int id) const
{
FPP sum = (FPP)0;
for(int i=0;i<this->getNbCol();i++)
sum += (*this)(id, i);
return sum;
}
template<typename FPP> template<typename FPP>
MatDense<FPP,Cpu> MatDense<FPP,Cpu>::get_block(faust_unsigned_int i, faust_unsigned_int j, faust_unsigned_int nrows, faust_unsigned_int ncols) MatDense<FPP,Cpu> MatDense<FPP,Cpu>::get_block(faust_unsigned_int i, faust_unsigned_int j, faust_unsigned_int nrows, faust_unsigned_int ncols)
{ {
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment