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 
176  (std::size_t m, std::size_t n, const scalar_t* D, std::size_t ld);
177 
194  (std::size_t m, std::size_t n, const DenseMatrix<scalar_t>& D,
195  std::size_t i, std::size_t j);
196 
198  DenseMatrix(const DenseMatrix<scalar_t>& D);
199 
201  DenseMatrix(DenseMatrix<scalar_t>&& D);
202 
204  virtual ~DenseMatrix();
205 
207  virtual DenseMatrix<scalar_t>& operator=(const DenseMatrix<scalar_t>& D);
208 
213  virtual DenseMatrix<scalar_t>& operator=(DenseMatrix<scalar_t>&& D);
214 
216  inline std::size_t rows() const { return rows_; }
217 
219  inline std::size_t cols() const { return cols_; }
220 
225  inline std::size_t ld() const { return ld_; }
226 
230  inline const scalar_t* data() const { return data_; }
231 
235  inline scalar_t* data() { return data_; }
236 
241  inline scalar_t* end() { return data_ + ld_ * cols_; }
242 
251  inline const scalar_t& operator()(std::size_t i, std::size_t j) const
252  { assert(i<=rows() && j<=cols()); return data_[i+ld_*j]; }
253 
262  inline const scalar_t* ptr(std::size_t i, std::size_t j) const
263  { assert(i<=rows() && j<=cols()); return data_+i+ld_*j; }
264 
273  inline scalar_t& operator()(std::size_t i, std::size_t j)
274  { assert(i<=rows() && j<=cols()); return data_[i+ld_*j]; }
275 
284  inline scalar_t* ptr(std::size_t i, std::size_t j)
285  { assert(i<=rows() && j<=cols()); return data_+i+ld_*j; }
286 
293  void print() const { print("A"); }
294 
308  void print(std::string name, bool all=false, int width=8) const;
309 
319  void print_to_file
320  (std::string name, std::string filename, int width=8) const;
321 
327  void random();
328 
333  void random
334  (random::RandomGeneratorBase<typename RealType<scalar_t>::
335  value_type>& rgen);
336 
338  void fill(scalar_t v);
339 
341  void zero();
342 
347  void eye();
348 
352  virtual void clear();
353 
362  void resize(std::size_t m, std::size_t n);
363 
372  void hconcat(const DenseMatrix<scalar_t>& b);
373 
385  void copy
386  (const DenseMatrix<scalar_t>& B, std::size_t i=0, std::size_t j=0);
387 
396  void copy(const scalar_t* B, std::size_t ldb);
397 
400 
404  void transpose(DenseMatrix<scalar_t>& X) const;
405 
416  void laswp(const std::vector<int>& P, bool fwd);
417 
428  void laswp(const int* P, bool fwd);
429 
442  void lapmr(const std::vector<int>& P, bool fwd);
443 
456  void lapmt(const std::vector<int>& P, bool fwd);
457 
468  void extract_rows
469  (const std::vector<std::size_t>& I, DenseMatrix<scalar_t>& B) const;
470 
481  (const std::vector<std::size_t>& I) const;
482 
493  void extract_cols
494  (const std::vector<std::size_t>& I, DenseMatrix<scalar_t>& B) const;
495 
506  (const std::vector<std::size_t>& I) const;
507 
519  (const std::vector<std::size_t>& I,
520  const std::vector<std::size_t>& J) const;
521 
533  (const std::vector<std::size_t>& I,
534  const DenseMatrix<scalar_t>& B, int depth);
535 
542  DenseMatrix<scalar_t>& add(const DenseMatrix<scalar_t>& B, int depth=0);
543 
551  DenseMatrix<scalar_t>& sub(const DenseMatrix<scalar_t>& B, int depth=0);
552 
558  DenseMatrix<scalar_t>& scale(scalar_t alpha, int depth=0);
559 
569  (scalar_t alpha, const DenseMatrix<scalar_t>& B, int depth=0);
570 
580  (scalar_t alpha, const DenseMatrix<scalar_t>& B, int depth=0);
581 
590  scale_rows(const std::vector<scalar_t>& D, int depth=0);
591 
593  scale_rows_real(const std::vector<real_t>& D, int depth=0);
594 
595 
604  DenseMatrix<scalar_t>& scale_rows(const scalar_t* D, int depth=0);
605 
606  DenseMatrix<scalar_t>& scale_rows_real(const real_t* D, int depth=0);
607 
616  (const std::vector<scalar_t>& D, int depth=0);
617 
622  real_t norm() const;
623 
627  real_t normF() const;
628 
632  real_t norm1() const;
633 
637  real_t normI() const;
638 
644  virtual std::size_t memory() const {
645  return sizeof(scalar_t) * rows() * cols();
646  }
647 
652  virtual std::size_t nonzeros() const {
653  return rows()*cols();
654  }
655 
668  std::vector<int> LU(int depth=0);
669 
679  int Cholesky(int depth=0);
680 
690  std::vector<int> LDLt(int depth=0);
691 
704  //std::vector<int> LDLt_rook(int depth=0);
705 
706 
720  (const DenseMatrix<scalar_t>& b,
721  const std::vector<int>& piv, int depth=0) const;
722 
734  void solve_LU_in_place
735  (DenseMatrix<scalar_t>& b, const std::vector<int>& piv, int depth=0) const;
736 
748  void solve_LU_in_place
749  (DenseMatrix<scalar_t>& b, const int* piv, int depth=0) const;
750 
763  (DenseMatrix<scalar_t>& b, const std::vector<int>& piv, int depth=0) const;
764 
779  // void solve_LDLt_rook_in_place
780  // (DenseMatrix<scalar_t>& b, const std::vector<int>& piv, int depth=0) const;
781 
794  void LQ
795  (DenseMatrix<scalar_t>& L, DenseMatrix<scalar_t>& Q, int depth) const;
796 
810  void orthogonalize(scalar_t& r_max, scalar_t& r_min, int depth);
811 
835  void ID_column
836  (DenseMatrix<scalar_t>& X, std::vector<int>& piv,
837  std::vector<std::size_t>& ind, real_t rel_tol,
838  real_t abs_tol, int max_rank, int depth);
839 
857  void ID_row
858  (DenseMatrix<scalar_t>& X, std::vector<int>& piv,
859  std::vector<std::size_t>& ind, real_t rel_tol, real_t abs_tol,
860  int max_rank, int depth) const;
861 
880  void low_rank
882  real_t rel_tol, real_t abs_tol, int max_rank, int depth) const;
883 
888  std::vector<scalar_t> singular_values() const;
889 
896  void shift(scalar_t sigma);
897 
898 
913  int syev(Jobz job, UpLo ul, std::vector<scalar_t>& lambda);
914 
915  private:
916  void ID_column_GEQP3
917  (DenseMatrix<scalar_t>& X, std::vector<int>& piv,
918  std::vector<std::size_t>& ind, real_t rel_tol,
919  real_t abs_tol, int max_rank, int depth);
920 
921  template<typename T> friend class DistributedMatrix;
922  };
923 
924 
945  template<typename scalar_t>
946  class DenseMatrixWrapper : public DenseMatrix<scalar_t> {
947  public:
951  DenseMatrixWrapper() : DenseMatrix<scalar_t>() {}
952 
964  (std::size_t m, std::size_t n, scalar_t* D, std::size_t ld) {
965  this->data_ = D; this->rows_ = m; this->cols_ = n;
966  this->ld_ = std::max(std::size_t(1), ld);
967  }
968 
983  (std::size_t m, std::size_t n, DenseMatrix<scalar_t>& D,
984  std::size_t i, std::size_t j)
985  : DenseMatrixWrapper<scalar_t>(m, n, &D(i, j), D.ld()) {
986  assert(i+m <= D.rows());
987  assert(j+n <= D.cols());
988  }
989 
995  virtual ~DenseMatrixWrapper() { this->data_ = nullptr; }
996 
1003  void clear() override {
1004  this->rows_ = 0; this->cols_ = 0;
1005  this->ld_ = 1; this->data_ = nullptr;
1006  }
1007 
1017  std::size_t memory() const override { return 0; }
1018 
1028  std::size_t nonzeros() const override { return 0; }
1029 
1034 
1040  DenseMatrixWrapper(const DenseMatrix<scalar_t>&) = delete;
1041 
1046 
1051 
1052  // /**
1053  // * Assignment operator. Shallow copy only. This only copies the
1054  // * wrapper object. Does not copy matrix elements.
1055  // *
1056  // * \param D matrix wrapper to copy from, this will be duplicated
1057  // */
1058  // DenseMatrixWrapper<scalar_t>&
1059  // operator=(const DenseMatrixWrapper<scalar_t>& D) {
1060  // this->data_ = D.data();
1061  // this->rows_ = D.rows();
1062  // this->cols_ = D.cols();
1063  // this->ld_ = D.ld();
1064  // return *this;
1065  // }
1066 
1075  this->data_ = D.data(); this->rows_ = D.rows();
1076  this->cols_ = D.cols(); this->ld_ = D.ld(); return *this; }
1077 
1087  operator=(const DenseMatrix<scalar_t>& a) override {
1088  assert(a.rows()==this->rows() && a.cols()==this->cols());
1089  for (std::size_t j=0; j<this->cols(); j++)
1090  for (std::size_t i=0; i<this->rows(); i++)
1091  this->operator()(i, j) = a(i, j);
1092  return *this;
1093  }
1094  };
1095 
1096 
1109  template<typename scalar_t>
1110  std::unique_ptr<const DenseMatrixWrapper<scalar_t>>
1112  (std::size_t m, std::size_t n, const scalar_t* D, std::size_t ld) {
1113  return std::unique_ptr<const DenseMatrixWrapper<scalar_t>>
1114  (new DenseMatrixWrapper<scalar_t>(m, n, const_cast<scalar_t*>(D), ld));
1115  }
1116 
1135  template<typename scalar_t>
1136  std::unique_ptr<const DenseMatrixWrapper<scalar_t>>
1138  (std::size_t m, std::size_t n, const DenseMatrix<scalar_t>& D,
1139  std::size_t i, std::size_t j) {
1140  return std::unique_ptr<const DenseMatrixWrapper<scalar_t>>
1142  (m, n, const_cast<DenseMatrix<scalar_t>&>(D), i, j));
1143  }
1144 
1145 
1162  template<typename scalar_t> void
1163  copy(std::size_t m, std::size_t n, const DenseMatrix<scalar_t>& a,
1164  std::size_t ia, std::size_t ja, DenseMatrix<scalar_t>& b,
1165  std::size_t ib, std::size_t jb) {
1166  for (std::size_t j=0; j<n; j++)
1167  for (std::size_t i=0; i<m; i++)
1168  b(ib+i, jb+j) = a(ia+i, ja+j);
1169  }
1170 
1181  template<typename scalar_t> void
1183  std::size_t ib, std::size_t jb) {
1184  copy(a.rows(), a.cols(), a, 0, 0, b, ib, jb);
1185  }
1186 
1195  template<typename scalar_t> void
1196  copy(const DenseMatrix<scalar_t>& a, scalar_t* b, std::size_t ldb) {
1197  for (std::size_t j=0; j<a.cols(); j++)
1198  for (std::size_t i=0; i<a.rows(); i++)
1199  b[i+j*ldb] = a(i, j);
1200  }
1201 
1202 
1211  template<typename scalar_t> DenseMatrix<scalar_t>
1213  assert(a.cols() == b.cols());
1214  DenseMatrix<scalar_t> tmp(a.rows()+b.rows(), a.cols());
1215  copy(a, tmp, 0, 0);
1216  copy(b, tmp, a.rows(), 0);
1217  return tmp;
1218  }
1219 
1228  template<typename scalar_t> DenseMatrix<scalar_t>
1230  assert(a.rows() == b.rows());
1231  DenseMatrix<scalar_t> tmp(a.rows(), a.cols()+b.cols());
1232  copy(a, tmp, 0, 0);
1233  copy(b, tmp, 0, a.cols());
1234  return tmp;
1235  }
1236 
1244  template<typename scalar_t> DenseMatrix<scalar_t>
1245  eye(std::size_t m, std::size_t n) {
1246  DenseMatrix<scalar_t> I(m, n);
1247  I.eye();
1248  return I;
1249  }
1250 
1251 
1252 
1253 
1270  template<typename scalar_t> void
1271  gemm(Trans ta, Trans tb, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1272  const DenseMatrix<scalar_t>& b, scalar_t beta,
1273  DenseMatrix<scalar_t>& c, int depth=0);
1274 
1275  template<typename scalar_t> void
1276  gemm(Trans ta, Trans tb, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1277  const scalar_t* b, int ldb, scalar_t beta,
1278  DenseMatrix<scalar_t>& c, int depth=0);
1279 
1280  template<typename scalar_t> void
1281  gemm(Trans ta, Trans tb, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1282  const DenseMatrix<scalar_t>& b, scalar_t beta,
1283  scalar_t* c, int ldc, int depth=0);
1284 
1294  template<typename scalar_t> void
1295  trmm(Side s, UpLo ul, Trans ta, Diag d, scalar_t alpha,
1296  const DenseMatrix<scalar_t>& a, DenseMatrix<scalar_t>& b,
1297  int depth=0);
1298 
1312  template<typename scalar_t> void
1313  trsm(Side s, UpLo ul, Trans ta, Diag d, scalar_t alpha,
1314  const DenseMatrix<scalar_t>& a, DenseMatrix<scalar_t>& b,
1315  int depth=0);
1316 
1325  template<typename scalar_t> void
1326  trsv(UpLo ul, Trans ta, Diag d, const DenseMatrix<scalar_t>& a,
1327  DenseMatrix<scalar_t>& b, int depth=0);
1328 
1337  template<typename scalar_t> void
1338  gemv(Trans ta, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1339  const DenseMatrix<scalar_t>& x, scalar_t beta,
1340  DenseMatrix<scalar_t>& y, int depth=0);
1341 
1350  template<typename scalar_t> void
1351  gemv(Trans ta, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1352  const scalar_t* x, int incx, scalar_t beta,
1353  DenseMatrix<scalar_t>& y, int depth=0);
1354 
1363  template<typename scalar_t> void
1364  gemv(Trans ta, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1365  const DenseMatrix<scalar_t>& x, scalar_t beta,
1366  scalar_t* y, int incy, int depth=0);
1367 
1376  template<typename scalar_t> void
1377  gemv(Trans ta, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1378  const scalar_t* x, int incx, scalar_t beta,
1379  scalar_t* y, int incy, int depth=0);
1380 
1381 
1383  template<typename scalar_t> long long int
1385  return (is_complex<scalar_t>() ? 4:1) *
1386  blas::getrf_flops(a.rows(), a.cols());
1387  }
1388 
1390  template<typename scalar_t> long long int
1392  return (is_complex<scalar_t>() ? 4:1) *
1393  blas::getrs_flops(b.rows(), b.cols());
1394  }
1395 
1397  template<typename scalar_t> long long int
1399  auto minrc = std::min(a.rows(), a.cols());
1400  return (is_complex<scalar_t>() ? 4:1) *
1401  (blas::gelqf_flops(a.rows(), a.cols()) +
1402  blas::xxglq_flops(a.cols(), a.cols(), minrc));
1403  }
1404 
1406  template<typename scalar_t> long long int
1407  ID_row_flops(const DenseMatrix<scalar_t>& a, int rank) {
1408  return (is_complex<scalar_t>() ? 4:1) *
1409  (blas::geqp3_flops(a.cols(), a.rows()) +
1410  blas::trsm_flops(rank, a.cols() - rank, scalar_t(1.), 'L'));
1411  }
1412 
1414  template<typename scalar_t> long long int
1415  trsm_flops(Side s, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1416  const DenseMatrix<scalar_t>& b) {
1417  return (is_complex<scalar_t>() ? 4:1) *
1418  blas::trsm_flops(b.rows(), b.cols(), alpha, char(s));
1419  }
1420 
1422  template<typename scalar_t> long long int
1423  gemm_flops(Trans ta, Trans tb, scalar_t alpha,
1424  const DenseMatrix<scalar_t>& a,
1425  const DenseMatrix<scalar_t>& b, scalar_t beta) {
1426  return (is_complex<scalar_t>() ? 4:1) *
1428  ((ta==Trans::N) ? a.rows() : a.cols(),
1429  (tb==Trans::N) ? b.cols() : b.rows(),
1430  (ta==Trans::N) ? a.cols() : a.rows(), alpha, beta);
1431  }
1432 
1434  template<typename scalar_t> long long int
1435  gemm_flops(Trans ta, Trans tb, scalar_t alpha,
1436  const DenseMatrix<scalar_t>& a, scalar_t beta,
1437  const DenseMatrix<scalar_t>& c) {
1438  return (is_complex<scalar_t>() ? 4:1) *
1440  (c.rows(), c.cols(), (ta==Trans::N) ? a.cols() : a.rows(), alpha, beta);
1441  }
1442 
1444  template<typename scalar_t> long long int
1446  auto minrc = std::min(a.rows(), a.cols());
1447  return (is_complex<scalar_t>() ? 4:1) *
1448  (blas::geqrf_flops(a.rows(), minrc) +
1449  blas::xxgqr_flops(a.rows(), minrc, minrc));
1450  }
1451 
1452 } // end namespace strumpack
1453 
1454 #endif // DENSE_MATRIX_HPP
strumpack::DenseMatrix::shift
void shift(scalar_t sigma)
strumpack::DenseMatrixWrapper::operator=
DenseMatrix< scalar_t > & operator=(const DenseMatrix< scalar_t > &a) override
Definition: DenseMatrix.hpp:1087
strumpack::ID_row_flops
long long int ID_row_flops(const DenseMatrix< scalar_t > &a, int rank)
Definition: DenseMatrix.hpp:1407
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:251
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:1435
strumpack::vconcat
DenseMatrix< scalar_t > vconcat(const DenseMatrix< scalar_t > &a, const DenseMatrix< scalar_t > &b)
Definition: DenseMatrix.hpp:1212
strumpack::DenseMatrix::cols
std::size_t cols() const
Definition: DenseMatrix.hpp:219
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:1415
strumpack::DenseMatrix::~DenseMatrix
virtual ~DenseMatrix()
strumpack::LQ_flops
long long int LQ_flops(const DenseMatrix< scalar_t > &a)
Definition: DenseMatrix.hpp:1398
strumpack::DenseMatrix::data
const scalar_t * data() const
Definition: DenseMatrix.hpp:230
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:1112
strumpack::Trans::C
@ C
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:995
strumpack::Jobz
Jobz
Definition: DenseMatrix.hpp:100
strumpack::DenseMatrix::nonzeros
virtual std::size_t nonzeros() const
Definition: DenseMatrix.hpp:652
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:1074
strumpack
Definition: StrumpackOptions.hpp:42
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:216
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:1423
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:946
strumpack::DenseMatrixWrapper::nonzeros
std::size_t nonzeros() const override
Definition: DenseMatrix.hpp:1028
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:1229
strumpack::DenseMatrix::operator=
virtual DenseMatrix< scalar_t > & operator=(const DenseMatrix< scalar_t > &D)
strumpack::DenseMatrix::ld
std::size_t ld() const
Definition: DenseMatrix.hpp:225
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:262
strumpack::DenseMatrix::ptr
scalar_t * ptr(std::size_t i, std::size_t j)
Definition: DenseMatrix.hpp:284
strumpack::DenseMatrix::scaled_add
DenseMatrix< scalar_t > & scaled_add(scalar_t alpha, const DenseMatrix< scalar_t > &B, int depth=0)
strumpack::copy
void copy(std::size_t m, std::size_t n, const DenseMatrix< scalar_t > &a, std::size_t ia, std::size_t ja, DenseMatrix< scalar_t > &b, std::size_t ib, std::size_t jb)
Definition: DenseMatrix.hpp:1163
strumpack::orthogonalize_flops
long long int orthogonalize_flops(const DenseMatrix< scalar_t > &a)
Definition: DenseMatrix.hpp:1445
strumpack::DenseMatrix::hconcat
void hconcat(const DenseMatrix< scalar_t > &b)
strumpack::UpLo
UpLo
Definition: DenseMatrix.hpp:82
strumpack::DenseMatrix::print
void print() const
Definition: DenseMatrix.hpp:293
strumpack::DenseMatrix::norm
real_t norm() const
strumpack::DenseMatrix::data
scalar_t * data()
Definition: DenseMatrix.hpp:235
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:241
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:644
strumpack::DenseMatrixWrapper::clear
void clear() override
Definition: DenseMatrix.hpp:1003
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:273
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:1245
strumpack::DenseMatrix::resize
void resize(std::size_t m, std::size_t n)
strumpack::DenseMatrix::LU
std::vector< int > LU(int depth=0)
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:1384
strumpack::DenseMatrix::DenseMatrix
DenseMatrix()
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::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:951
strumpack::solve_flops
long long int solve_flops(const DenseMatrix< scalar_t > &b)
Definition: DenseMatrix.hpp:1391
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:1017