DenseMatrix.hpp
Go to the documentation of this file.
1 /*
2  * STRUMPACK -- STRUctured Matrices PACKage, Copyright (c) 2014, The
3  * Regents of the University of California, through Lawrence Berkeley
4  * National Laboratory (subject to receipt of any required approvals
5  * from the U.S. Dept. of Energy). All rights reserved.
6  *
7  * If you have questions about your rights to use or distribute this
8  * software, please contact Berkeley Lab's Technology Transfer
9  * Department at TTD@lbl.gov.
10  *
11  * NOTICE. This software is owned by the U.S. Department of Energy. As
12  * such, the U.S. Government has been granted for itself and others
13  * acting on its behalf a paid-up, nonexclusive, irrevocable,
14  * worldwide license in the Software to reproduce, prepare derivative
15  * works, and perform publicly and display publicly. Beginning five
16  * (5) years after the date permission to assert copyright is obtained
17  * from the U.S. Department of Energy, and subject to any subsequent
18  * five (5) year renewals, the U.S. Government is granted for itself
19  * and others acting on its behalf a paid-up, nonexclusive,
20  * irrevocable, worldwide license in the Software to reproduce,
21  * prepare derivative works, distribute copies to the public, perform
22  * publicly and display publicly, and to permit others to do so.
23  *
24  * Developers: Pieter Ghysels, Francois-Henry Rouet, Xiaoye S. Li.
25  * (Lawrence Berkeley National Lab, Computational Research
26  * Division).
27  *
28  */
34 #ifndef DENSE_MATRIX_HPP
35 #define DENSE_MATRIX_HPP
36 
37 #include <string>
38 #include <vector>
39 
40 #include "misc/RandomWrapper.hpp"
41 #include "BLASLAPACKWrapper.hpp"
42 
43 
44 namespace strumpack {
45 
50  enum class Trans : char {
51  N='N',
52  C='C',
53  T='T'
54  };
55 
56  inline Trans c2T(char op) {
57  switch (op) {
58  case 'n': case 'N': return Trans::N;
59  case 't': case 'T': return Trans::T;
60  case 'c': case 'C': return Trans::C;
61  default:
62  std::cerr << "ERROR: char " << op << " not recognized,"
63  << " should be one of n/N, t/T or c/C" << std::endl;
64  return Trans::N;
65  }
66  }
67 
68 
73  enum class Side : char {
74  L='L',
75  R='R'
76  };
77 
82  enum class UpLo : char {
83  U='U',
84  L='L'
85  };
86 
91  enum class Diag : char {
92  U='U',
93  N='N'
94  };
95 
100  enum class Jobz : char {
101  N='N',
102  V='V'
103  };
104 
105 
137  template<typename scalar_t> class DenseMatrix {
138  using real_t = typename RealType<scalar_t>::value_type;
139 
140  protected:
141  scalar_t* data_ = nullptr;
142  std::size_t rows_ = 0;
143  std::size_t cols_ = 0;
144  std::size_t ld_ = 1;
145 
146  public:
147 
152  DenseMatrix();
153 
161  DenseMatrix(std::size_t m, std::size_t n);
162 
175  DenseMatrix(std::size_t m, std::size_t n,
176  const scalar_t* D, std::size_t ld);
177 
193  DenseMatrix(std::size_t m, std::size_t n, const DenseMatrix<scalar_t>& D,
194  std::size_t i, std::size_t j);
195 
197  DenseMatrix(const DenseMatrix<scalar_t>& D);
198 
200  DenseMatrix(DenseMatrix<scalar_t>&& D);
201 
203  virtual ~DenseMatrix();
204 
206  virtual DenseMatrix<scalar_t>& operator=(const DenseMatrix<scalar_t>& D);
207 
212  virtual DenseMatrix<scalar_t>& operator=(DenseMatrix<scalar_t>&& D);
213 
215  inline std::size_t rows() const { return rows_; }
216 
218  inline std::size_t cols() const { return cols_; }
219 
224  inline std::size_t ld() const { return ld_; }
225 
229  inline const scalar_t* data() const { return data_; }
230 
234  inline scalar_t* data() { return data_; }
235 
240  inline scalar_t* end() { return data_ + ld_ * cols_; }
241 
246  inline const scalar_t* end() const {
247  if (rows_ == 0 || cols_ == 0) return data_;
248  return data_ + ld_ * cols_;
249  }
250 
259  inline const scalar_t& operator()(std::size_t i, std::size_t j) const
260  { assert(i<=rows() && j<=cols()); return data_[i+ld_*j]; }
261 
270  inline const scalar_t* ptr(std::size_t i, std::size_t j) const
271  { assert(i<=rows() && j<=cols()); return data_+i+ld_*j; }
272 
281  inline scalar_t& operator()(std::size_t i, std::size_t j)
282  { assert(i<=rows() && j<=cols()); return data_[i+ld_*j]; }
283 
292  inline scalar_t* ptr(std::size_t i, std::size_t j)
293  { assert(i<=rows() && j<=cols()); return data_+i+ld_*j; }
294 
301  void print() const { print("A"); }
302 
316  void print(std::string name, bool all=false, int width=8) const;
317 
327  void print_to_file(std::string name,
328  std::string filename, int width=8) const;
329 
335  void random();
336 
341  void random(random::RandomGeneratorBase<typename RealType<scalar_t>::
342  value_type>& rgen);
343 
345  void fill(scalar_t v);
346 
348  void zero();
349 
354  void eye();
355 
359  virtual void clear();
360 
369  void resize(std::size_t m, std::size_t n);
370 
379  void hconcat(const DenseMatrix<scalar_t>& b);
380 
392  void copy(const DenseMatrix<scalar_t>& B,
393  std::size_t i=0, std::size_t j=0);
394 
403  void copy(const scalar_t* B, std::size_t ldb);
404 
407 
411  void transpose(DenseMatrix<scalar_t>& X) const;
412 
423  void laswp(const std::vector<int>& P, bool fwd);
424 
435  void laswp(const int* P, bool fwd);
436 
449  void lapmr(const std::vector<int>& P, bool fwd);
450 
463  void lapmt(const std::vector<int>& P, bool fwd);
464 
475  void extract_rows(const std::vector<std::size_t>& I,
476  DenseMatrix<scalar_t>& B) const;
477 
488  extract_rows(const std::vector<std::size_t>& I) const;
489 
500  void extract_cols(const std::vector<std::size_t>& I,
501  DenseMatrix<scalar_t>& B) const;
502 
513  extract_cols(const std::vector<std::size_t>& I) const;
514 
526  extract(const std::vector<std::size_t>& I,
527  const std::vector<std::size_t>& J) const;
528 
540  scatter_rows_add(const std::vector<std::size_t>& I,
541  const DenseMatrix<scalar_t>& B, int depth);
542 
549  DenseMatrix<scalar_t>& add(const DenseMatrix<scalar_t>& B, int depth=0);
550 
558  DenseMatrix<scalar_t>& sub(const DenseMatrix<scalar_t>& B, int depth=0);
559 
565  DenseMatrix<scalar_t>& scale(scalar_t alpha, int depth=0);
566 
576  scaled_add(scalar_t alpha, const DenseMatrix<scalar_t>& B, int depth=0);
577 
587  (scalar_t alpha, const DenseMatrix<scalar_t>& B, int depth=0);
588 
597  scale_rows(const std::vector<scalar_t>& D, int depth=0);
598 
600  scale_rows_real(const std::vector<real_t>& D, int depth=0);
601 
602 
611  DenseMatrix<scalar_t>& scale_rows(const scalar_t* D, int depth=0);
612 
613  DenseMatrix<scalar_t>& scale_rows_real(const real_t* D, int depth=0);
614 
623  (const std::vector<scalar_t>& D, int depth=0);
624 
629  real_t norm() const;
630 
634  real_t normF() const;
635 
639  real_t norm1() const;
640 
644  real_t normI() const;
645 
651  virtual std::size_t memory() const {
652  return sizeof(scalar_t) * rows() * cols();
653  }
654 
659  virtual std::size_t nonzeros() const {
660  return rows()*cols();
661  }
662 
677  int LU(std::vector<int>& piv, int depth=0);
678 
692  std::vector<int> LU(int depth=0);
693 
703  int Cholesky(int depth=0);
704 
714  std::vector<int> LDLt(int depth=0);
715 
728  //std::vector<int> LDLt_rook(int depth=0);
729 
730 
744  (const DenseMatrix<scalar_t>& b,
745  const std::vector<int>& piv, int depth=0) const;
746 
758  void solve_LU_in_place
759  (DenseMatrix<scalar_t>& b, const std::vector<int>& piv, int depth=0) const;
760 
772  void solve_LU_in_place
773  (DenseMatrix<scalar_t>& b, const int* piv, int depth=0) const;
774 
787  (DenseMatrix<scalar_t>& b, const std::vector<int>& piv, int depth=0) const;
788 
803  // void solve_LDLt_rook_in_place
804  // (DenseMatrix<scalar_t>& b, const std::vector<int>& piv, int depth=0) const;
805 
818  void LQ
819  (DenseMatrix<scalar_t>& L, DenseMatrix<scalar_t>& Q, int depth) const;
820 
834  void orthogonalize(scalar_t& r_max, scalar_t& r_min, int depth);
835 
859  void ID_column
860  (DenseMatrix<scalar_t>& X, std::vector<int>& piv,
861  std::vector<std::size_t>& ind, real_t rel_tol,
862  real_t abs_tol, int max_rank, int depth);
863 
881  void ID_row
882  (DenseMatrix<scalar_t>& X, std::vector<int>& piv,
883  std::vector<std::size_t>& ind, real_t rel_tol, real_t abs_tol,
884  int max_rank, int depth) const;
885 
904  void low_rank
906  real_t rel_tol, real_t abs_tol, int max_rank, int depth) const;
907 
912  std::vector<scalar_t> singular_values() const;
913 
920  void shift(scalar_t sigma);
921 
922 
937  int syev(Jobz job, UpLo ul, std::vector<scalar_t>& lambda);
938 
945  void write(const std::string& fname) const;
946 
953  static DenseMatrix<scalar_t> read(const std::string& fname);
954 
955  private:
956  void ID_column_GEQP3
957  (DenseMatrix<scalar_t>& X, std::vector<int>& piv,
958  std::vector<std::size_t>& ind, real_t rel_tol,
959  real_t abs_tol, int max_rank, int depth);
960 
961  template<typename T> friend class DistributedMatrix;
962 
963  template<typename T> friend std::ofstream&
964  operator<<(std::ofstream& os, const DenseMatrix<T>& D);
965  template<typename T> friend std::ifstream&
966  operator>>(std::ifstream& is, DenseMatrix<T>& D);
967  };
968 
969 
990  template<typename scalar_t>
991  class DenseMatrixWrapper : public DenseMatrix<scalar_t> {
992  public:
996  DenseMatrixWrapper() : DenseMatrix<scalar_t>() {}
997 
1008  DenseMatrixWrapper(std::size_t m, std::size_t n,
1009  scalar_t* D, std::size_t ld) {
1010  this->data_ = D; this->rows_ = m; this->cols_ = n;
1011  this->ld_ = std::max(std::size_t(1), ld);
1012  }
1013 
1027  DenseMatrixWrapper(std::size_t m, std::size_t n, DenseMatrix<scalar_t>& D,
1028  std::size_t i, std::size_t j)
1029  : DenseMatrixWrapper<scalar_t>(m, n, &D(i, j), D.ld()) {
1030  assert(i+m <= D.rows());
1031  assert(j+n <= D.cols());
1032  }
1033 
1039  virtual ~DenseMatrixWrapper() { this->data_ = nullptr; }
1040 
1047  void clear() override {
1048  this->rows_ = 0; this->cols_ = 0;
1049  this->ld_ = 1; this->data_ = nullptr;
1050  }
1051 
1061  std::size_t memory() const override { return 0; }
1062 
1072  std::size_t nonzeros() const override { return 0; }
1073 
1078 
1084  DenseMatrixWrapper(const DenseMatrix<scalar_t>&) = delete;
1085 
1090 
1095 
1096  // /**
1097  // * Assignment operator. Shallow copy only. This only copies the
1098  // * wrapper object. Does not copy matrix elements.
1099  // *
1100  // * \param D matrix wrapper to copy from, this will be duplicated
1101  // */
1102  // DenseMatrixWrapper<scalar_t>&
1103  // operator=(const DenseMatrixWrapper<scalar_t>& D) {
1104  // this->data_ = D.data();
1105  // this->rows_ = D.rows();
1106  // this->cols_ = D.cols();
1107  // this->ld_ = D.ld();
1108  // return *this;
1109  // }
1110 
1119  this->data_ = D.data(); this->rows_ = D.rows();
1120  this->cols_ = D.cols(); this->ld_ = D.ld(); return *this; }
1121 
1131  operator=(const DenseMatrix<scalar_t>& a) override {
1132  assert(a.rows()==this->rows() && a.cols()==this->cols());
1133  for (std::size_t j=0; j<this->cols(); j++)
1134  for (std::size_t i=0; i<this->rows(); i++)
1135  this->operator()(i, j) = a(i, j);
1136  return *this;
1137  }
1138  };
1139 
1140 
1153  template<typename scalar_t>
1154  std::unique_ptr<const DenseMatrixWrapper<scalar_t>>
1156  (std::size_t m, std::size_t n, const scalar_t* D, std::size_t ld) {
1157  return std::unique_ptr<const DenseMatrixWrapper<scalar_t>>
1158  (new DenseMatrixWrapper<scalar_t>(m, n, const_cast<scalar_t*>(D), ld));
1159  }
1160 
1179  template<typename scalar_t>
1180  std::unique_ptr<const DenseMatrixWrapper<scalar_t>>
1182  (std::size_t m, std::size_t n, const DenseMatrix<scalar_t>& D,
1183  std::size_t i, std::size_t j) {
1184  return std::unique_ptr<const DenseMatrixWrapper<scalar_t>>
1186  (m, n, const_cast<DenseMatrix<scalar_t>&>(D), i, j));
1187  }
1188 
1189 
1206  template<typename scalar_from_t, typename scalar_to_t> void
1207  copy(std::size_t m, std::size_t n, const DenseMatrix<scalar_from_t>& a,
1208  std::size_t ia, std::size_t ja, DenseMatrix<scalar_to_t>& b,
1209  std::size_t ib, std::size_t jb) {
1210  for (std::size_t j=0; j<n; j++)
1211  for (std::size_t i=0; i<m; i++)
1212  b(ib+i, jb+j) = static_cast<scalar_to_t>(a(ia+i, ja+j));
1213  }
1214 
1225  template<typename scalar_from_t, typename scalar_to_t> void
1227  std::size_t ib=0, std::size_t jb=0) {
1228  copy(a.rows(), a.cols(), a, 0, 0, b, ib, jb);
1229  }
1230 
1239  template<typename scalar_t> void
1240  copy(const DenseMatrix<scalar_t>& a, scalar_t* b, std::size_t ldb) {
1241  for (std::size_t j=0; j<a.cols(); j++)
1242  for (std::size_t i=0; i<a.rows(); i++)
1243  b[i+j*ldb] = a(i, j);
1244  }
1245 
1246 
1255  template<typename scalar_t> DenseMatrix<scalar_t>
1257  assert(a.cols() == b.cols());
1258  DenseMatrix<scalar_t> tmp(a.rows()+b.rows(), a.cols());
1259  copy(a, tmp, 0, 0);
1260  copy(b, tmp, a.rows(), 0);
1261  return tmp;
1262  }
1263 
1272  template<typename scalar_t> DenseMatrix<scalar_t>
1274  assert(a.rows() == b.rows());
1275  DenseMatrix<scalar_t> tmp(a.rows(), a.cols()+b.cols());
1276  copy(a, tmp, 0, 0);
1277  copy(b, tmp, 0, a.cols());
1278  return tmp;
1279  }
1280 
1288  template<typename scalar_t> DenseMatrix<scalar_t>
1289  eye(std::size_t m, std::size_t n) {
1290  DenseMatrix<scalar_t> I(m, n);
1291  I.eye();
1292  return I;
1293  }
1294 
1295 
1296 
1297 
1314  template<typename scalar_t> void
1315  gemm(Trans ta, Trans tb, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1316  const DenseMatrix<scalar_t>& b, scalar_t beta,
1317  DenseMatrix<scalar_t>& c, int depth=0);
1318 
1319  template<typename scalar_t> void
1320  gemm(Trans ta, Trans tb, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1321  const scalar_t* b, int ldb, scalar_t beta,
1322  DenseMatrix<scalar_t>& c, int depth=0);
1323 
1324  template<typename scalar_t> void
1325  gemm(Trans ta, Trans tb, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1326  const DenseMatrix<scalar_t>& b, scalar_t beta,
1327  scalar_t* c, int ldc, int depth=0);
1328 
1338  template<typename scalar_t> void
1339  trmm(Side s, UpLo ul, Trans ta, Diag d, scalar_t alpha,
1340  const DenseMatrix<scalar_t>& a, DenseMatrix<scalar_t>& b,
1341  int depth=0);
1342 
1356  template<typename scalar_t> void
1357  trsm(Side s, UpLo ul, Trans ta, Diag d, scalar_t alpha,
1358  const DenseMatrix<scalar_t>& a, DenseMatrix<scalar_t>& b,
1359  int depth=0);
1360 
1369  template<typename scalar_t> void
1370  trsv(UpLo ul, Trans ta, Diag d, const DenseMatrix<scalar_t>& a,
1371  DenseMatrix<scalar_t>& b, int depth=0);
1372 
1381  template<typename scalar_t> void
1382  gemv(Trans ta, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1383  const DenseMatrix<scalar_t>& x, scalar_t beta,
1384  DenseMatrix<scalar_t>& y, int depth=0);
1385 
1394  template<typename scalar_t> void
1395  gemv(Trans ta, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1396  const scalar_t* x, int incx, scalar_t beta,
1397  DenseMatrix<scalar_t>& y, int depth=0);
1398 
1407  template<typename scalar_t> void
1408  gemv(Trans ta, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1409  const DenseMatrix<scalar_t>& x, scalar_t beta,
1410  scalar_t* y, int incy, int depth=0);
1411 
1420  template<typename scalar_t> void
1421  gemv(Trans ta, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1422  const scalar_t* x, int incx, scalar_t beta,
1423  scalar_t* y, int incy, int depth=0);
1424 
1425 
1427  template<typename scalar_t> long long int
1429  return (is_complex<scalar_t>() ? 4:1) *
1430  blas::getrf_flops(a.rows(), a.cols());
1431  }
1432 
1434  template<typename scalar_t> long long int
1436  return (is_complex<scalar_t>() ? 4:1) *
1437  blas::getrs_flops(b.rows(), b.cols());
1438  }
1439 
1441  template<typename scalar_t> long long int
1443  auto minrc = std::min(a.rows(), a.cols());
1444  return (is_complex<scalar_t>() ? 4:1) *
1445  (blas::gelqf_flops(a.rows(), a.cols()) +
1446  blas::xxglq_flops(a.cols(), a.cols(), minrc));
1447  }
1448 
1450  template<typename scalar_t> long long int
1451  ID_row_flops(const DenseMatrix<scalar_t>& a, int rank) {
1452  return (is_complex<scalar_t>() ? 4:1) *
1453  (blas::geqp3_flops(a.cols(), a.rows()) +
1454  blas::trsm_flops(rank, a.cols() - rank, scalar_t(1.), 'L'));
1455  }
1456 
1458  template<typename scalar_t> long long int
1459  trsm_flops(Side s, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1460  const DenseMatrix<scalar_t>& b) {
1461  return (is_complex<scalar_t>() ? 4:1) *
1462  blas::trsm_flops(b.rows(), b.cols(), alpha, char(s));
1463  }
1464 
1466  template<typename scalar_t> long long int
1467  gemm_flops(Trans ta, Trans tb, scalar_t alpha,
1468  const DenseMatrix<scalar_t>& a,
1469  const DenseMatrix<scalar_t>& b, scalar_t beta) {
1470  return (is_complex<scalar_t>() ? 4:1) *
1472  ((ta==Trans::N) ? a.rows() : a.cols(),
1473  (tb==Trans::N) ? b.cols() : b.rows(),
1474  (ta==Trans::N) ? a.cols() : a.rows(), alpha, beta);
1475  }
1476 
1478  template<typename scalar_t> long long int
1479  gemm_flops(Trans ta, Trans tb, scalar_t alpha,
1480  const DenseMatrix<scalar_t>& a, scalar_t beta,
1481  const DenseMatrix<scalar_t>& c) {
1482  return (is_complex<scalar_t>() ? 4:1) *
1484  (c.rows(), c.cols(), (ta==Trans::N) ? a.cols() : a.rows(), alpha, beta);
1485  }
1486 
1488  template<typename scalar_t> long long int
1490  auto minrc = std::min(a.rows(), a.cols());
1491  return (is_complex<scalar_t>() ? 4:1) *
1492  (blas::geqrf_flops(a.rows(), minrc) +
1493  blas::xxgqr_flops(a.rows(), minrc, minrc));
1494  }
1495 
1505  template<typename scalar_t,typename cast_t>
1506  DenseMatrix<cast_t> cast_matrix(const DenseMatrix<scalar_t>& mat);
1507 
1508 } // end namespace strumpack
1509 
1510 #endif // DENSE_MATRIX_HPP
strumpack::DenseMatrix::shift
void shift(scalar_t sigma)
strumpack::DenseMatrix::LU
int LU(std::vector< int > &piv, int depth=0)
strumpack::DenseMatrixWrapper::operator=
DenseMatrix< scalar_t > & operator=(const DenseMatrix< scalar_t > &a) override
Definition: DenseMatrix.hpp:1131
strumpack::ID_row_flops
long long int ID_row_flops(const DenseMatrix< scalar_t > &a, int rank)
Definition: DenseMatrix.hpp:1451
strumpack::DenseMatrix::scale_and_add
DenseMatrix< scalar_t > & scale_and_add(scalar_t alpha, const DenseMatrix< scalar_t > &B, int depth=0)
strumpack::DenseMatrix::extract_cols
void extract_cols(const std::vector< std::size_t > &I, DenseMatrix< scalar_t > &B) const
strumpack::DenseMatrix::operator()
const scalar_t & operator()(std::size_t i, std::size_t j) const
Definition: DenseMatrix.hpp:259
strumpack::UpLo::U
@ U
strumpack::gemm_flops
long long int gemm_flops(Trans ta, Trans tb, scalar_t alpha, const DenseMatrix< scalar_t > &a, scalar_t beta, const DenseMatrix< scalar_t > &c)
Definition: DenseMatrix.hpp:1479
strumpack::vconcat
DenseMatrix< scalar_t > vconcat(const DenseMatrix< scalar_t > &a, const DenseMatrix< scalar_t > &b)
Definition: DenseMatrix.hpp:1256
strumpack::DenseMatrix::cols
std::size_t cols() const
Definition: DenseMatrix.hpp:218
strumpack::DenseMatrix::lapmt
void lapmt(const std::vector< int > &P, bool fwd)
strumpack::Trans::N
@ N
strumpack::trsm_flops
long long int trsm_flops(Side s, scalar_t alpha, const DenseMatrix< scalar_t > &a, const DenseMatrix< scalar_t > &b)
Definition: DenseMatrix.hpp:1459
strumpack::DenseMatrix::~DenseMatrix
virtual ~DenseMatrix()
strumpack::LQ_flops
long long int LQ_flops(const DenseMatrix< scalar_t > &a)
Definition: DenseMatrix.hpp:1442
strumpack::DenseMatrix::data
const scalar_t * data() const
Definition: DenseMatrix.hpp:229
strumpack::DenseMatrix::sub
DenseMatrix< scalar_t > & sub(const DenseMatrix< scalar_t > &B, int depth=0)
strumpack::ConstDenseMatrixWrapperPtr
std::unique_ptr< const DenseMatrixWrapper< scalar_t > > ConstDenseMatrixWrapperPtr(std::size_t m, std::size_t n, const scalar_t *D, std::size_t ld)
Definition: DenseMatrix.hpp:1156
strumpack::Trans::C
@ C
strumpack::copy
void copy(std::size_t m, std::size_t n, const DenseMatrix< scalar_from_t > &a, std::size_t ia, std::size_t ja, DenseMatrix< scalar_to_t > &b, std::size_t ib, std::size_t jb)
Definition: DenseMatrix.hpp:1207
strumpack::DenseMatrix::extract_rows
void extract_rows(const std::vector< std::size_t > &I, DenseMatrix< scalar_t > &B) const
strumpack::Side::R
@ R
strumpack::DenseMatrixWrapper::~DenseMatrixWrapper
virtual ~DenseMatrixWrapper()
Definition: DenseMatrix.hpp:1039
strumpack::Jobz
Jobz
Definition: DenseMatrix.hpp:100
strumpack::DenseMatrix::nonzeros
virtual std::size_t nonzeros() const
Definition: DenseMatrix.hpp:659
strumpack::DenseMatrix::lapmr
void lapmr(const std::vector< int > &P, bool fwd)
strumpack::random::RandomGeneratorBase
class to wrap the C++11 random number generator/distribution
Definition: RandomWrapper.hpp:105
strumpack::DenseMatrixWrapper::operator=
DenseMatrixWrapper< scalar_t > & operator=(DenseMatrixWrapper< scalar_t > &&D)
Definition: DenseMatrix.hpp:1118
strumpack
Definition: StrumpackOptions.hpp:42
strumpack::operator<<
std::ostream & operator<<(std::ostream &os, const BLACSGrid *g)
Definition: BLACSGrid.hpp:351
strumpack::trsv
void trsv(UpLo ul, Trans ta, Diag d, const DenseMatrix< scalar_t > &a, DenseMatrix< scalar_t > &b, int depth=0)
strumpack::DistributedMatrix
Definition: CompressedSparseMatrix.hpp:58
strumpack::DenseMatrix::rows
std::size_t rows() const
Definition: DenseMatrix.hpp:215
strumpack::gemm_flops
long long int gemm_flops(Trans ta, Trans tb, scalar_t alpha, const DenseMatrix< scalar_t > &a, const DenseMatrix< scalar_t > &b, scalar_t beta)
Definition: DenseMatrix.hpp:1467
strumpack::DenseMatrix::normF
real_t normF() const
strumpack::Diag::N
@ N
strumpack::DenseMatrix::Cholesky
int Cholesky(int depth=0)
strumpack::DenseMatrix::LDLt
std::vector< int > LDLt(int depth=0)
strumpack::DenseMatrixWrapper
Like DenseMatrix, this class represents a matrix, stored in column major format, to allow direct use ...
Definition: DenseMatrix.hpp:991
strumpack::DenseMatrixWrapper::nonzeros
std::size_t nonzeros() const override
Definition: DenseMatrix.hpp:1072
strumpack::trsm
void trsm(Side s, UpLo ul, Trans ta, Diag d, scalar_t alpha, const DenseMatrix< scalar_t > &a, DenseMatrix< scalar_t > &b, int depth=0)
strumpack::DenseMatrix::copy
void copy(const DenseMatrix< scalar_t > &B, std::size_t i=0, std::size_t j=0)
strumpack::DenseMatrix< scalar_t >
strumpack::DenseMatrix::LQ
void LQ(DenseMatrix< scalar_t > &L, DenseMatrix< scalar_t > &Q, int depth) const
strumpack::DenseMatrix::transpose
DenseMatrix< scalar_t > transpose() const
strumpack::Jobz::N
@ N
strumpack::hconcat
DenseMatrix< scalar_t > hconcat(const DenseMatrix< scalar_t > &a, const DenseMatrix< scalar_t > &b)
Definition: DenseMatrix.hpp:1273
strumpack::DenseMatrix::operator=
virtual DenseMatrix< scalar_t > & operator=(const DenseMatrix< scalar_t > &D)
strumpack::DenseMatrix::ld
std::size_t ld() const
Definition: DenseMatrix.hpp:224
strumpack::DenseMatrix::ID_column
void ID_column(DenseMatrix< scalar_t > &X, std::vector< int > &piv, std::vector< std::size_t > &ind, real_t rel_tol, real_t abs_tol, int max_rank, int depth)
strumpack::DenseMatrix::ptr
const scalar_t * ptr(std::size_t i, std::size_t j) const
Definition: DenseMatrix.hpp:270
strumpack::DenseMatrix::ptr
scalar_t * ptr(std::size_t i, std::size_t j)
Definition: DenseMatrix.hpp:292
strumpack::DenseMatrix::scaled_add
DenseMatrix< scalar_t > & scaled_add(scalar_t alpha, const DenseMatrix< scalar_t > &B, int depth=0)
strumpack::orthogonalize_flops
long long int orthogonalize_flops(const DenseMatrix< scalar_t > &a)
Definition: DenseMatrix.hpp:1489
strumpack::DenseMatrix::hconcat
void hconcat(const DenseMatrix< scalar_t > &b)
strumpack::UpLo
UpLo
Definition: DenseMatrix.hpp:82
strumpack::DenseMatrix::read
static DenseMatrix< scalar_t > read(const std::string &fname)
strumpack::DenseMatrixWrapper::DenseMatrixWrapper
DenseMatrixWrapper(std::size_t m, std::size_t n, scalar_t *D, std::size_t ld)
Definition: DenseMatrix.hpp:1008
strumpack::DenseMatrix::print
void print() const
Definition: DenseMatrix.hpp:301
strumpack::DenseMatrix::end
const scalar_t * end() const
Definition: DenseMatrix.hpp:246
strumpack::DenseMatrix::norm
real_t norm() const
strumpack::DenseMatrix::data
scalar_t * data()
Definition: DenseMatrix.hpp:234
strumpack::DenseMatrix::add
DenseMatrix< scalar_t > & add(const DenseMatrix< scalar_t > &B, int depth=0)
strumpack::DenseMatrix::laswp
void laswp(const std::vector< int > &P, bool fwd)
strumpack::DenseMatrix::extract
DenseMatrix< scalar_t > extract(const std::vector< std::size_t > &I, const std::vector< std::size_t > &J) const
strumpack::DenseMatrix::end
scalar_t * end()
Definition: DenseMatrix.hpp:240
strumpack::DenseMatrix::scatter_rows_add
DenseMatrix< scalar_t > & scatter_rows_add(const std::vector< std::size_t > &I, const DenseMatrix< scalar_t > &B, int depth)
strumpack::Side::L
@ L
strumpack::DenseMatrix::memory
virtual std::size_t memory() const
Definition: DenseMatrix.hpp:651
strumpack::DenseMatrixWrapper::clear
void clear() override
Definition: DenseMatrix.hpp:1047
strumpack::DenseMatrix::solve_LU_in_place
void solve_LU_in_place(DenseMatrix< scalar_t > &b, const std::vector< int > &piv, int depth=0) const
strumpack::DenseMatrix::random
void random()
strumpack::DenseMatrix::orthogonalize
void orthogonalize(scalar_t &r_max, scalar_t &r_min, int depth)
strumpack::DenseMatrix::operator()
scalar_t & operator()(std::size_t i, std::size_t j)
Definition: DenseMatrix.hpp:281
strumpack::DenseMatrix::normI
real_t normI() const
strumpack::gemv
void gemv(Trans ta, scalar_t alpha, const DenseMatrix< scalar_t > &a, const DenseMatrix< scalar_t > &x, scalar_t beta, DenseMatrix< scalar_t > &y, int depth=0)
strumpack::eye
DenseMatrix< scalar_t > eye(std::size_t m, std::size_t n)
Definition: DenseMatrix.hpp:1289
strumpack::DenseMatrix::resize
void resize(std::size_t m, std::size_t n)
strumpack::DenseMatrix::write
void write(const std::string &fname) const
strumpack::trmm
void trmm(Side s, UpLo ul, Trans ta, Diag d, scalar_t alpha, const DenseMatrix< scalar_t > &a, DenseMatrix< scalar_t > &b, int depth=0)
strumpack::DenseMatrix::syev
int syev(Jobz job, UpLo ul, std::vector< scalar_t > &lambda)
strumpack::DenseMatrix::ID_row
void ID_row(DenseMatrix< scalar_t > &X, std::vector< int > &piv, std::vector< std::size_t > &ind, real_t rel_tol, real_t abs_tol, int max_rank, int depth) const
strumpack::gemm
void gemm(Trans ta, Trans tb, scalar_t alpha, const DenseMatrix< scalar_t > &a, const DenseMatrix< scalar_t > &b, scalar_t beta, DenseMatrix< scalar_t > &c, int depth=0)
strumpack::DenseMatrix::clear
virtual void clear()
strumpack::DenseMatrix::fill
void fill(scalar_t v)
strumpack::DenseMatrix::zero
void zero()
strumpack::UpLo::L
@ L
strumpack::DenseMatrix::scale
DenseMatrix< scalar_t > & scale(scalar_t alpha, int depth=0)
strumpack::Side
Side
Definition: DenseMatrix.hpp:73
strumpack::DenseMatrix::singular_values
std::vector< scalar_t > singular_values() const
strumpack::LU_flops
long long int LU_flops(const DenseMatrix< scalar_t > &a)
Definition: DenseMatrix.hpp:1428
strumpack::DenseMatrix::DenseMatrix
DenseMatrix()
strumpack::cast_matrix
CSRMatrix< cast_t, integer_t > cast_matrix(const CSRMatrix< scalar_t, integer_t > &mat)
strumpack::Diag::U
@ U
strumpack::DenseMatrix::scale_rows
DenseMatrix< scalar_t > & scale_rows(const std::vector< scalar_t > &D, int depth=0)
strumpack::DenseMatrix::div_rows
DenseMatrix< scalar_t > & div_rows(const std::vector< scalar_t > &D, int depth=0)
strumpack::DenseMatrixWrapper::DenseMatrixWrapper
DenseMatrixWrapper(std::size_t m, std::size_t n, DenseMatrix< scalar_t > &D, std::size_t i, std::size_t j)
Definition: DenseMatrix.hpp:1027
strumpack::Diag
Diag
Definition: DenseMatrix.hpp:91
strumpack::DenseMatrix::norm1
real_t norm1() const
strumpack::DenseMatrix::low_rank
void low_rank(DenseMatrix< scalar_t > &U, DenseMatrix< scalar_t > &V, real_t rel_tol, real_t abs_tol, int max_rank, int depth) const
strumpack::DenseMatrixWrapper::DenseMatrixWrapper
DenseMatrixWrapper()
Definition: DenseMatrix.hpp:996
strumpack::solve_flops
long long int solve_flops(const DenseMatrix< scalar_t > &b)
Definition: DenseMatrix.hpp:1435
strumpack::Jobz::V
@ V
strumpack::Trans::T
@ T
strumpack::DenseMatrix::eye
void eye()
strumpack::DenseMatrix::solve
DenseMatrix< scalar_t > solve(const DenseMatrix< scalar_t > &b, const std::vector< int > &piv, int depth=0) const
strumpack::DenseMatrix::solve_LDLt_in_place
void solve_LDLt_in_place(DenseMatrix< scalar_t > &b, const std::vector< int > &piv, int depth=0) const
strumpack::Trans
Trans
Definition: DenseMatrix.hpp:50
strumpack::DenseMatrix::print_to_file
void print_to_file(std::string name, std::string filename, int width=8) const
strumpack::DenseMatrixWrapper::memory
std::size_t memory() const override
Definition: DenseMatrix.hpp:1061