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#include <functional>
40
41#include "misc/RandomWrapper.hpp"
42#include "BLASLAPACKWrapper.hpp"
43
44
45namespace strumpack {
46
51 enum class Trans : char {
52 N='N',
53 C='C',
54 T='T'
55 };
56
57 inline Trans c2T(char op) {
58 switch (op) {
59 case 'n': case 'N': return Trans::N;
60 case 't': case 'T': return Trans::T;
61 case 'c': case 'C': return Trans::C;
62 default:
63 std::cerr << "ERROR: char " << op << " not recognized,"
64 << " should be one of n/N, t/T or c/C" << std::endl;
65 return Trans::N;
66 }
67 }
68
69
74 enum class Side : char {
75 L='L',
76 R='R'
77 };
78
83 enum class UpLo : char {
84 U='U',
85 L='L'
86 };
87
92 enum class Diag : char {
93 U='U',
94 N='N'
95 };
96
101 enum class Jobz : char {
102 N='N',
103 V='V'
104 };
105
106
138 template<typename scalar_t> class DenseMatrix {
139 using real_t = typename RealType<scalar_t>::value_type;
140
141 protected:
142 scalar_t* data_ = nullptr;
143 std::size_t rows_ = 0;
144 std::size_t cols_ = 0;
145 std::size_t ld_ = 1;
146
147 public:
148
154
162 DenseMatrix(std::size_t m, std::size_t n);
163
172 DenseMatrix(std::size_t m, std::size_t n,
173 const std::function<scalar_t(std::size_t,std::size_t)>& A);
174
187 DenseMatrix(std::size_t m, std::size_t n,
188 const scalar_t* D, std::size_t ld);
189
205 DenseMatrix(std::size_t m, std::size_t n, const DenseMatrix<scalar_t>& D,
206 std::size_t i, std::size_t j);
207
210
213
215 virtual ~DenseMatrix();
216
219
225
227 inline std::size_t rows() const { return rows_; }
228
230 inline std::size_t cols() const { return cols_; }
231
236 inline std::size_t ld() const { return ld_; }
237
241 inline const scalar_t* data() const { return data_; }
242
246 inline scalar_t* data() { return data_; }
247
252 inline scalar_t* end() { return data_ + ld_ * cols_; }
253
258 inline const scalar_t* end() const {
259 if (rows_ == 0 || cols_ == 0) return data_;
260 return data_ + ld_ * cols_;
261 }
262
271 inline const scalar_t& operator()(std::size_t i, std::size_t j) const
272 { assert(i<=rows() && j<=cols()); return data_[i+ld_*j]; }
273
282 inline const scalar_t* ptr(std::size_t i, std::size_t j) const
283 { assert(i<=rows() && j<=cols()); return data_+i+ld_*j; }
284
293 inline scalar_t& operator()(std::size_t i, std::size_t j)
294 { assert(i<=rows() && j<=cols()); return data_[i+ld_*j]; }
295
304 inline scalar_t* ptr(std::size_t i, std::size_t j)
305 { assert(i<=rows() && j<=cols()); return data_+i+ld_*j; }
306
313 void print() const { print("A"); }
314
328 void print(std::string name, bool all=false, int width=8) const;
329
339 void print_to_file(std::string name,
340 std::string filename, int width=8) const;
341
347 void random();
348
353 void random(random::RandomGeneratorBase<typename RealType<scalar_t>::
354 value_type>& rgen);
355
361 void fill(scalar_t v);
362
369 void fill(const std::function<scalar_t(std::size_t,std::size_t)>& A);
370
372 void zero();
373
378 void eye();
379
383 virtual void clear();
384
393 void resize(std::size_t m, std::size_t n);
394
404
417 std::size_t i=0, std::size_t j=0);
418
427 void copy(const scalar_t* B, std::size_t ldb);
428
431
436
447 void laswp(const std::vector<int>& P, bool fwd);
448
459 void laswp(const int* P, bool fwd);
460
473 void lapmr(const std::vector<int>& P, bool fwd);
474
487 void lapmt(const std::vector<int>& P, bool fwd);
488
499 void extract_rows(const std::vector<std::size_t>& I,
500 DenseMatrix<scalar_t>& B) const;
501
512 extract_rows(const std::vector<std::size_t>& I) const;
513
524 void extract_cols(const std::vector<std::size_t>& I,
525 DenseMatrix<scalar_t>& B) const;
526
537 extract_cols(const std::vector<std::size_t>& I) const;
538
550 extract(const std::vector<std::size_t>& I,
551 const std::vector<std::size_t>& J) const;
552
564 scatter_rows_add(const std::vector<std::size_t>& I,
565 const DenseMatrix<scalar_t>& B, int depth);
566
574
583
589 DenseMatrix<scalar_t>& scale(scalar_t alpha, int depth=0);
590
600 scaled_add(scalar_t alpha, const DenseMatrix<scalar_t>& B, int depth=0);
601
611 (scalar_t alpha, const DenseMatrix<scalar_t>& B, int depth=0);
612
621 scale_rows(const std::vector<scalar_t>& D, int depth=0);
622
624 scale_rows_real(const std::vector<real_t>& D, int depth=0);
625
626
635 DenseMatrix<scalar_t>& scale_rows(const scalar_t* D, int depth=0);
636
637 DenseMatrix<scalar_t>& scale_rows_real(const real_t* D, int depth=0);
638
647 (const std::vector<scalar_t>& D, int depth=0);
648
653 real_t norm() const;
654
658 real_t normF() const;
659
663 real_t norm1() const;
664
668 real_t normI() const;
669
675 virtual std::size_t memory() const {
676 return sizeof(scalar_t) * rows() * cols();
677 }
678
683 virtual std::size_t nonzeros() const {
684 return rows()*cols();
685 }
686
701 int LU(std::vector<int>& piv, int depth=0);
702
716 std::vector<int> LU(int depth=0);
717
727 int Cholesky(int depth=0);
728
738 std::vector<int> LDLt(int depth=0);
739
752 //std::vector<int> LDLt_rook(int depth=0);
753
754
768 (const DenseMatrix<scalar_t>& b,
769 const std::vector<int>& piv, int depth=0) const;
770
783 (DenseMatrix<scalar_t>& b, const std::vector<int>& piv, int depth=0) const;
784
797 (DenseMatrix<scalar_t>& b, const int* piv, int depth=0) const;
798
811 (DenseMatrix<scalar_t>& b, const std::vector<int>& piv, int depth=0) const;
812
827 // void solve_LDLt_rook_in_place
828 // (DenseMatrix<scalar_t>& b, const std::vector<int>& piv, int depth=0) const;
829
842 void LQ
843 (DenseMatrix<scalar_t>& L, DenseMatrix<scalar_t>& Q, int depth) const;
844
858 void orthogonalize(scalar_t& r_max, scalar_t& r_min, int depth);
859
884 (DenseMatrix<scalar_t>& X, std::vector<int>& piv,
885 std::vector<std::size_t>& ind, real_t rel_tol,
886 real_t abs_tol, int max_rank, int depth);
887
906 (DenseMatrix<scalar_t>& X, std::vector<int>& piv,
907 std::vector<std::size_t>& ind, real_t rel_tol, real_t abs_tol,
908 int max_rank, int depth) const;
909
930 real_t rel_tol, real_t abs_tol, int max_rank, int depth) const;
931
936 std::vector<scalar_t> singular_values() const;
937
944 void shift(scalar_t sigma);
945
946
961 int syev(Jobz job, UpLo ul, std::vector<scalar_t>& lambda);
962
969 void write(const std::string& fname) const;
970
977 static DenseMatrix<scalar_t> read(const std::string& fname);
978
979 private:
980 void ID_column_GEQP3
981 (DenseMatrix<scalar_t>& X, std::vector<int>& piv,
982 std::vector<std::size_t>& ind, real_t rel_tol,
983 real_t abs_tol, int max_rank, int depth);
984
985 template<typename T> friend class DistributedMatrix;
986
987 template<typename T> friend std::ofstream&
988 operator<<(std::ofstream& os, const DenseMatrix<T>& D);
989 template<typename T> friend std::ifstream&
990 operator>>(std::ifstream& is, DenseMatrix<T>& D);
991 };
992
993
1014 template<typename scalar_t>
1015 class DenseMatrixWrapper : public DenseMatrix<scalar_t> {
1016 public:
1021
1032 DenseMatrixWrapper(std::size_t m, std::size_t n,
1033 scalar_t* D, std::size_t ld) {
1034 this->data_ = D; this->rows_ = m; this->cols_ = n;
1035 this->ld_ = std::max(std::size_t(1), ld);
1036 }
1037
1051 DenseMatrixWrapper(std::size_t m, std::size_t n, DenseMatrix<scalar_t>& D,
1052 std::size_t i, std::size_t j)
1053 : DenseMatrixWrapper<scalar_t>(m, n, &D(i, j), D.ld()) {
1054 assert(i+m <= D.rows());
1055 assert(j+n <= D.cols());
1056 }
1057
1063 virtual ~DenseMatrixWrapper() { this->data_ = nullptr; }
1064
1071 void clear() override {
1072 this->rows_ = 0; this->cols_ = 0;
1073 this->ld_ = 1; this->data_ = nullptr;
1074 }
1075
1085 std::size_t memory() const override { return 0; }
1086
1096 std::size_t nonzeros() const override { return 0; }
1097
1102
1109
1114
1119
1120 // /**
1121 // * Assignment operator. Shallow copy only. This only copies the
1122 // * wrapper object. Does not copy matrix elements.
1123 // *
1124 // * \param D matrix wrapper to copy from, this will be duplicated
1125 // */
1126 // DenseMatrixWrapper<scalar_t>&
1127 // operator=(const DenseMatrixWrapper<scalar_t>& D) {
1128 // this->data_ = D.data();
1129 // this->rows_ = D.rows();
1130 // this->cols_ = D.cols();
1131 // this->ld_ = D.ld();
1132 // return *this;
1133 // }
1134
1143 this->data_ = D.data(); this->rows_ = D.rows();
1144 this->cols_ = D.cols(); this->ld_ = D.ld(); return *this; }
1145
1156 assert(a.rows()==this->rows() && a.cols()==this->cols());
1157 for (std::size_t j=0; j<this->cols(); j++)
1158 for (std::size_t i=0; i<this->rows(); i++)
1159 this->operator()(i, j) = a(i, j);
1160 return *this;
1161 }
1162 };
1163
1164
1177 template<typename scalar_t>
1178 std::unique_ptr<const DenseMatrixWrapper<scalar_t>>
1180 (std::size_t m, std::size_t n, const scalar_t* D, std::size_t ld) {
1181 return std::unique_ptr<const DenseMatrixWrapper<scalar_t>>
1182 (new DenseMatrixWrapper<scalar_t>(m, n, const_cast<scalar_t*>(D), ld));
1183 }
1184
1203 template<typename scalar_t>
1204 std::unique_ptr<const DenseMatrixWrapper<scalar_t>>
1206 (std::size_t m, std::size_t n, const DenseMatrix<scalar_t>& D,
1207 std::size_t i, std::size_t j) {
1208 return std::unique_ptr<const DenseMatrixWrapper<scalar_t>>
1210 (m, n, const_cast<DenseMatrix<scalar_t>&>(D), i, j));
1211 }
1212
1213
1230 template<typename scalar_from_t, typename scalar_to_t> void
1231 copy(std::size_t m, std::size_t n, const DenseMatrix<scalar_from_t>& a,
1232 std::size_t ia, std::size_t ja, DenseMatrix<scalar_to_t>& b,
1233 std::size_t ib, std::size_t jb) {
1234 for (std::size_t j=0; j<n; j++)
1235 for (std::size_t i=0; i<m; i++)
1236 b(ib+i, jb+j) = static_cast<scalar_to_t>(a(ia+i, ja+j));
1237 }
1238
1249 template<typename scalar_from_t, typename scalar_to_t> void
1251 std::size_t ib=0, std::size_t jb=0) {
1252 copy(a.rows(), a.cols(), a, 0, 0, b, ib, jb);
1253 }
1254
1263 template<typename scalar_t> void
1264 copy(const DenseMatrix<scalar_t>& a, scalar_t* b, std::size_t ldb) {
1265 for (std::size_t j=0; j<a.cols(); j++)
1266 for (std::size_t i=0; i<a.rows(); i++)
1267 b[i+j*ldb] = a(i, j);
1268 }
1269
1270
1279 template<typename scalar_t> DenseMatrix<scalar_t>
1281 assert(a.cols() == b.cols());
1282 DenseMatrix<scalar_t> tmp(a.rows()+b.rows(), a.cols());
1283 copy(a, tmp, 0, 0);
1284 copy(b, tmp, a.rows(), 0);
1285 return tmp;
1286 }
1287
1296 template<typename scalar_t> DenseMatrix<scalar_t>
1298 assert(a.rows() == b.rows());
1299 DenseMatrix<scalar_t> tmp(a.rows(), a.cols()+b.cols());
1300 copy(a, tmp, 0, 0);
1301 copy(b, tmp, 0, a.cols());
1302 return tmp;
1303 }
1304
1312 template<typename scalar_t> DenseMatrix<scalar_t>
1313 eye(std::size_t m, std::size_t n) {
1314 DenseMatrix<scalar_t> I(m, n);
1315 I.eye();
1316 return I;
1317 }
1318
1319
1320
1321
1338 template<typename scalar_t> void
1339 gemm(Trans ta, Trans tb, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1340 const DenseMatrix<scalar_t>& b, scalar_t beta,
1341 DenseMatrix<scalar_t>& c, int depth=0);
1342
1343 template<typename scalar_t> void
1344 gemm(Trans ta, Trans tb, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1345 const scalar_t* b, int ldb, scalar_t beta,
1346 DenseMatrix<scalar_t>& c, int depth=0);
1347
1348 template<typename scalar_t> void
1349 gemm(Trans ta, Trans tb, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1350 const DenseMatrix<scalar_t>& b, scalar_t beta,
1351 scalar_t* c, int ldc, int depth=0);
1352
1362 template<typename scalar_t> void
1363 trmm(Side s, UpLo ul, Trans ta, Diag d, scalar_t alpha,
1365 int depth=0);
1366
1380 template<typename scalar_t> void
1381 trsm(Side s, UpLo ul, Trans ta, Diag d, scalar_t alpha,
1383 int depth=0);
1384
1393 template<typename scalar_t> void
1395 DenseMatrix<scalar_t>& b, int depth=0);
1396
1405 template<typename scalar_t> void
1406 gemv(Trans ta, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1407 const DenseMatrix<scalar_t>& x, scalar_t beta,
1408 DenseMatrix<scalar_t>& y, int depth=0);
1409
1418 template<typename scalar_t> void
1419 gemv(Trans ta, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1420 const scalar_t* x, int incx, scalar_t beta,
1421 DenseMatrix<scalar_t>& y, int depth=0);
1422
1431 template<typename scalar_t> void
1432 gemv(Trans ta, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1433 const DenseMatrix<scalar_t>& x, scalar_t beta,
1434 scalar_t* y, int incy, int depth=0);
1435
1444 template<typename scalar_t> void
1445 gemv(Trans ta, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1446 const scalar_t* x, int incx, scalar_t beta,
1447 scalar_t* y, int incy, int depth=0);
1448
1449
1451 template<typename scalar_t> long long int
1453 return (is_complex<scalar_t>() ? 4:1) *
1454 blas::getrf_flops(a.rows(), a.cols());
1455 }
1456
1458 template<typename scalar_t> long long int
1460 return (is_complex<scalar_t>() ? 4:1) *
1461 blas::getrs_flops(b.rows(), b.cols());
1462 }
1463
1465 template<typename scalar_t> long long int
1467 auto minrc = std::min(a.rows(), a.cols());
1468 return (is_complex<scalar_t>() ? 4:1) *
1469 (blas::gelqf_flops(a.rows(), a.cols()) +
1470 blas::xxglq_flops(a.cols(), a.cols(), minrc));
1471 }
1472
1474 template<typename scalar_t> long long int
1476 return (is_complex<scalar_t>() ? 4:1) *
1477 (blas::geqp3_flops(a.cols(), a.rows()) +
1478 blas::trsm_flops(rank, a.cols() - rank, scalar_t(1.), 'L'));
1479 }
1480
1482 template<typename scalar_t> long long int
1483 trsm_flops(Side s, scalar_t alpha, const DenseMatrix<scalar_t>& a,
1484 const DenseMatrix<scalar_t>& b) {
1485 return (is_complex<scalar_t>() ? 4:1) *
1486 blas::trsm_flops(b.rows(), b.cols(), alpha, char(s));
1487 }
1488
1490 template<typename scalar_t> long long int
1491 gemm_flops(Trans ta, Trans tb, scalar_t alpha,
1492 const DenseMatrix<scalar_t>& a,
1493 const DenseMatrix<scalar_t>& b, scalar_t beta) {
1494 return (is_complex<scalar_t>() ? 4:1) *
1496 ((ta==Trans::N) ? a.rows() : a.cols(),
1497 (tb==Trans::N) ? b.cols() : b.rows(),
1498 (ta==Trans::N) ? a.cols() : a.rows(), alpha, beta);
1499 }
1500
1502 template<typename scalar_t> long long int
1503 gemm_flops(Trans ta, Trans tb, scalar_t alpha,
1504 const DenseMatrix<scalar_t>& a, scalar_t beta,
1505 const DenseMatrix<scalar_t>& c) {
1506 return (is_complex<scalar_t>() ? 4:1) *
1508 (c.rows(), c.cols(), (ta==Trans::N) ? a.cols() : a.rows(), alpha, beta);
1509 }
1510
1512 template<typename scalar_t> long long int
1514 auto minrc = std::min(a.rows(), a.cols());
1515 return (is_complex<scalar_t>() ? 4:1) *
1516 (blas::geqrf_flops(a.rows(), minrc) +
1517 blas::xxgqr_flops(a.rows(), minrc, minrc));
1518 }
1519
1529 template<typename scalar_t,typename cast_t>
1531
1532} // end namespace strumpack
1533
1534#endif // DENSE_MATRIX_HPP
Like DenseMatrix, this class represents a matrix, stored in column major format, to allow direct use ...
Definition: DenseMatrix.hpp:1015
DenseMatrix< scalar_t > & operator=(const DenseMatrix< scalar_t > &a) override
Definition: DenseMatrix.hpp:1155
virtual ~DenseMatrixWrapper()
Definition: DenseMatrix.hpp:1063
DenseMatrixWrapper(const DenseMatrixWrapper< scalar_t > &)=default
DenseMatrixWrapper()
Definition: DenseMatrix.hpp:1020
DenseMatrixWrapper(const DenseMatrix< scalar_t > &)=delete
DenseMatrixWrapper(std::size_t m, std::size_t n, scalar_t *D, std::size_t ld)
Definition: DenseMatrix.hpp:1032
void clear() override
Definition: DenseMatrix.hpp:1071
std::size_t nonzeros() const override
Definition: DenseMatrix.hpp:1096
DenseMatrixWrapper(std::size_t m, std::size_t n, DenseMatrix< scalar_t > &D, std::size_t i, std::size_t j)
Definition: DenseMatrix.hpp:1051
std::size_t memory() const override
Definition: DenseMatrix.hpp:1085
DenseMatrixWrapper(DenseMatrixWrapper< scalar_t > &&)=default
DenseMatrixWrapper(DenseMatrix< scalar_t > &&)=delete
DenseMatrixWrapper< scalar_t > & operator=(DenseMatrixWrapper< scalar_t > &&D)
Definition: DenseMatrix.hpp:1142
This class represents a matrix, stored in column major format, to allow direct use of BLAS/LAPACK rou...
Definition: DenseMatrix.hpp:138
real_t norm() const
void lapmt(const std::vector< int > &P, bool fwd)
DenseMatrix< scalar_t > solve(const DenseMatrix< scalar_t > &b, const std::vector< int > &piv, int depth=0) const
static DenseMatrix< scalar_t > read(const std::string &fname)
real_t norm1() const
DenseMatrix(DenseMatrix< scalar_t > &&D)
void copy(const scalar_t *B, std::size_t ldb)
void solve_LDLt_in_place(DenseMatrix< scalar_t > &b, const std::vector< int > &piv, int depth=0) const
std::size_t cols() const
Definition: DenseMatrix.hpp:230
DenseMatrix< scalar_t > & scale_and_add(scalar_t alpha, const DenseMatrix< scalar_t > &B, int depth=0)
void laswp(const int *P, bool fwd)
scalar_t * ptr(std::size_t i, std::size_t j)
Definition: DenseMatrix.hpp:304
void fill(const std::function< scalar_t(std::size_t, std::size_t)> &A)
virtual std::size_t nonzeros() const
Definition: DenseMatrix.hpp:683
void random(random::RandomGeneratorBase< typename RealType< scalar_t >::value_type > &rgen)
void copy(const DenseMatrix< scalar_t > &B, std::size_t i=0, std::size_t j=0)
DenseMatrix< scalar_t > & sub(const DenseMatrix< scalar_t > &B, int depth=0)
DenseMatrix< scalar_t > & scale_rows(const scalar_t *D, int depth=0)
DenseMatrix(std::size_t m, std::size_t n, const scalar_t *D, std::size_t ld)
std::vector< scalar_t > singular_values() const
std::size_t rows() const
Definition: DenseMatrix.hpp:227
std::size_t ld() const
Definition: DenseMatrix.hpp:236
DenseMatrix(std::size_t m, std::size_t n, const DenseMatrix< scalar_t > &D, std::size_t i, std::size_t j)
virtual DenseMatrix< scalar_t > & operator=(DenseMatrix< scalar_t > &&D)
void print(std::string name, bool all=false, int width=8) const
void laswp(const std::vector< int > &P, bool fwd)
void fill(scalar_t v)
void solve_LU_in_place(DenseMatrix< scalar_t > &b, const std::vector< int > &piv, int depth=0) const
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)
const scalar_t * ptr(std::size_t i, std::size_t j) const
Definition: DenseMatrix.hpp:282
scalar_t * data()
Definition: DenseMatrix.hpp:246
void lapmr(const std::vector< int > &P, bool fwd)
DenseMatrix(const DenseMatrix< scalar_t > &D)
void LQ(DenseMatrix< scalar_t > &L, DenseMatrix< scalar_t > &Q, int depth) const
virtual DenseMatrix< scalar_t > & operator=(const DenseMatrix< scalar_t > &D)
void print_to_file(std::string name, std::string filename, int width=8) const
scalar_t * end()
Definition: DenseMatrix.hpp:252
void extract_cols(const std::vector< std::size_t > &I, DenseMatrix< scalar_t > &B) const
DenseMatrix< scalar_t > extract_cols(const std::vector< std::size_t > &I) const
DenseMatrix< scalar_t > extract_rows(const std::vector< std::size_t > &I) const
real_t normI() const
void hconcat(const DenseMatrix< scalar_t > &b)
std::vector< int > LDLt(int depth=0)
real_t normF() const
virtual std::size_t memory() const
Definition: DenseMatrix.hpp:675
void solve_LU_in_place(DenseMatrix< scalar_t > &b, const int *piv, int depth=0) const
int syev(Jobz job, UpLo ul, std::vector< scalar_t > &lambda)
DenseMatrix< scalar_t > & add(const DenseMatrix< scalar_t > &B, int depth=0)
DenseMatrix< scalar_t > & div_rows(const std::vector< scalar_t > &D, int depth=0)
DenseMatrix< scalar_t > & scale(scalar_t alpha, int depth=0)
DenseMatrix< scalar_t > & scatter_rows_add(const std::vector< std::size_t > &I, const DenseMatrix< scalar_t > &B, int depth)
void transpose(DenseMatrix< scalar_t > &X) const
scalar_t & operator()(std::size_t i, std::size_t j)
Definition: DenseMatrix.hpp:293
DenseMatrix< scalar_t > extract(const std::vector< std::size_t > &I, const std::vector< std::size_t > &J) const
DenseMatrix(std::size_t m, std::size_t n)
void extract_rows(const std::vector< std::size_t > &I, DenseMatrix< scalar_t > &B) const
void orthogonalize(scalar_t &r_max, scalar_t &r_min, int depth)
int Cholesky(int depth=0)
DenseMatrix< scalar_t > & scaled_add(scalar_t alpha, const DenseMatrix< scalar_t > &B, int depth=0)
DenseMatrix< scalar_t > & scale_rows(const std::vector< scalar_t > &D, int depth=0)
int LU(std::vector< int > &piv, int depth=0)
void write(const std::string &fname) const
virtual void clear()
void print() const
Definition: DenseMatrix.hpp:313
std::vector< int > LU(int depth=0)
const scalar_t * end() const
Definition: DenseMatrix.hpp:258
DenseMatrix< scalar_t > transpose() const
void resize(std::size_t m, std::size_t n)
const scalar_t & operator()(std::size_t i, std::size_t j) const
Definition: DenseMatrix.hpp:271
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
const scalar_t * data() const
Definition: DenseMatrix.hpp:241
DenseMatrix(std::size_t m, std::size_t n, const std::function< scalar_t(std::size_t, std::size_t)> &A)
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
void shift(scalar_t sigma)
2D block cyclicly distributed matrix, as used by ScaLAPACK.
Definition: DistributedMatrix.hpp:84
class to wrap the C++11 random number generator/distribution
Definition: RandomWrapper.hpp:105
Definition: StrumpackOptions.hpp:43
long long int LU_flops(const DenseMatrix< scalar_t > &a)
Definition: DenseMatrix.hpp:1452
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:1491
DenseMatrix< scalar_t > eye(std::size_t m, std::size_t n)
Definition: DenseMatrix.hpp:1313
UpLo
Definition: DenseMatrix.hpp:83
CSRMatrix< cast_t, integer_t > cast_matrix(const CSRMatrix< scalar_t, integer_t > &mat)
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:1180
DenseMatrix< scalar_t > vconcat(const DenseMatrix< scalar_t > &a, const DenseMatrix< scalar_t > &b)
Definition: DenseMatrix.hpp:1280
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)
void trsv(UpLo ul, Trans ta, Diag d, const DenseMatrix< scalar_t > &a, DenseMatrix< scalar_t > &b, int depth=0)
DenseMatrix< scalar_t > hconcat(const DenseMatrix< scalar_t > &a, const DenseMatrix< scalar_t > &b)
Definition: DenseMatrix.hpp:1297
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)
Trans
Definition: DenseMatrix.hpp:51
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)
long long int orthogonalize_flops(const DenseMatrix< scalar_t > &a)
Definition: DenseMatrix.hpp:1513
long long int trsm_flops(Side s, scalar_t alpha, const DenseMatrix< scalar_t > &a, const DenseMatrix< scalar_t > &b)
Definition: DenseMatrix.hpp:1483
Side
Definition: DenseMatrix.hpp:74
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:1231
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)
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:1503
long long int ID_row_flops(const DenseMatrix< scalar_t > &a, int rank)
Definition: DenseMatrix.hpp:1475
long long int LQ_flops(const DenseMatrix< scalar_t > &a)
Definition: DenseMatrix.hpp:1466
long long int solve_flops(const DenseMatrix< scalar_t > &b)
Definition: DenseMatrix.hpp:1459
Diag
Definition: DenseMatrix.hpp:92
Jobz
Definition: DenseMatrix.hpp:101