10 static_assert(std::is_same_v<T, double> ||
11 std::is_same_v<T, std::complex<double>>,
12 "Determinant only works for double");
14 assert(rows() == cols() &&
"Determinant only defined for square matrix");
18 auto gsl_view = LU.as_gsl_view();
19 gsl_permutation *permutn = gsl_permutation_alloc(rows());
20 if constexpr (std::is_same_v<T, double>) {
21 gsl_linalg_LU_decomp(&gsl_view.matrix, permutn, &sLU);
22 gsl_permutation_free(permutn);
23 return gsl_linalg_LU_det(&gsl_view.matrix, sLU);
24 }
else if constexpr (std::is_same_v<T, std::complex<double>>) {
25 gsl_linalg_complex_LU_decomp(&gsl_view.matrix, permutn, &sLU);
26 gsl_permutation_free(permutn);
27 const auto gsl_cmplx = gsl_linalg_complex_LU_det(&gsl_view.matrix, sLU);
29 return {GSL_REAL(gsl_cmplx), GSL_IMAG(gsl_cmplx)};
39 std::is_same_v<T, double> || std::is_same_v<T, std::complex<double>>,
40 "invert only works for Matrix<double> or Matrix<complex<double>>");
42 assert(rows() == cols() &&
"Inverse only defined for square matrix");
49 auto LU_gsl = LU.as_gsl_view();
50 auto iverse_gsl = this->as_gsl_view();
51 gsl_permutation *permutn = gsl_permutation_alloc(m_rows);
52 if constexpr (std::is_same_v<T, double>) {
53 gsl_linalg_LU_decomp(&LU_gsl.matrix, permutn, &sLU);
54 gsl_linalg_LU_invert(&LU_gsl.matrix, permutn, &iverse_gsl.matrix);
55 }
else if constexpr (std::is_same_v<T, std::complex<double>>) {
56 gsl_linalg_complex_LU_decomp(&LU_gsl.matrix, permutn, &sLU);
57 gsl_linalg_complex_LU_invert(&LU_gsl.matrix, permutn, &iverse_gsl.matrix);
59 gsl_permutation_free(permutn);
67 if constexpr (std::is_same_v<T, double>) {
69 const auto this_gsl = as_gsl_view();
70 gsl_matrix_transpose_memcpy(&Tr_gsl.matrix, &this_gsl.matrix);
71 }
else if constexpr (std::is_same_v<T, float>) {
73 const auto this_gsl = as_gsl_view();
74 gsl_matrix_float_transpose_memcpy(&Tr_gsl.matrix, &this_gsl.matrix);
75 }
else if constexpr (std::is_same_v<T, std::complex<double>>) {
77 const auto this_gsl = as_gsl_view();
78 gsl_matrix_complex_transpose_memcpy(&Tr_gsl.matrix, &this_gsl.matrix);
79 }
else if constexpr (std::is_same_v<T, std::complex<float>>) {
81 const auto this_gsl = as_gsl_view();
82 gsl_matrix_complex_float_transpose_memcpy(&Tr_gsl.matrix, &this_gsl.matrix);
85 for (
auto i = 0ul; i < Tr.
rows(); ++i) {
86 for (
auto j = 0ul; j < Tr.
cols(); ++j) {
87 Tr[i][j] = (*this)[j][i];
98 assert(m_rows == m_cols &&
"Can only call make_identity() for square matrix");
99 for (
auto i = 0ul; i < m_rows; ++i) {
100 for (
auto j = 0ul; j < m_cols; ++j) {
101 at(i, j) = i == j ? T(1) : T(0);
109 for (std::size_t i = 0; i < size(); ++i) {
118 static_assert(is_complex_v<T>,
"conj() only available for complex Matrix");
119 std::vector<T> conj_data;
120 conj_data.reserve(m_data.size());
121 for (std::size_t i = 0; i < m_data.size(); ++i) {
122 conj_data.push_back(std::conj(m_data[i]));
124 return Matrix<T>{m_rows, m_cols, std::move(conj_data)};
129 static_assert(is_complex_v<T>,
"conj() only available for complex Matrix");
130 for (std::size_t i = 0; i < m_data.size(); ++i) {
131 m_data[i] = std::conj(m_data[i]);
138 static_assert(is_complex_v<T>,
"real() only available for complex Matrix");
139 std::vector<typename T::value_type> real_data;
140 real_data.reserve(m_data.size());
141 for (std::size_t i = 0; i < m_data.size(); ++i) {
142 real_data.push_back(std::real(m_data[i]));
149 static_assert(is_complex_v<T>,
"imag() only available for complex Matrix");
150 std::vector<typename T::value_type> imag_data;
151 imag_data.reserve(m_data.size());
152 for (std::size_t i = 0; i < m_data.size(); ++i) {
153 imag_data.push_back(std::imag(m_data[i]));
160 static_assert(!is_complex_v<T>,
"complex() only available for real Matrix");
162 std::vector<std::complex<T>> new_data;
163 new_data.reserve(m_data.size());
164 for (std::size_t i = 0; i < m_data.size(); ++i) {
165 new_data.push_back(m_data[i]);
173 assert(rows() == rhs.
rows() && cols() == rhs.
cols() &&
174 "Matrices must have same dimensions for addition");
176 this->m_data += rhs.m_data;
181 assert(rows() == rhs.
rows() && cols() == rhs.
cols() &&
182 "Matrices must have same dimensions for subtraction");
184 this->m_data -= rhs.m_data;
205 assert(m_rows == m_cols &&
"Can only call M+a for square matrix");
206 for (
auto i = 0ul; i < m_rows; ++i) {
215 assert(m_rows == m_cols &&
"Can only call M-a for square matrix");
216 for (
auto i = 0ul; i < m_rows; ++i) {
225 assert(rows() == a.
rows() && cols() == a.
cols() &&
226 "Matrices must have same dimensions for mult_elements_by");
227 for (
auto i = 0ul; i < m_data.size(); ++i) {
228 m_data[i] *= a.m_data[i];
238 "Matrices a and b must have correct dimension for multiplication");
241 GEMM(a, b, &product);
258 const int A_rows =
static_cast<int>(trans_A ? a.
cols() : a.
rows());
259 const int A_cols =
static_cast<int>(trans_A ? a.
rows() : a.
cols());
260 const int B_rows =
static_cast<int>(trans_B ? b.
cols() : b.
rows());
261 const int B_cols =
static_cast<int>(trans_B ? b.
rows() : b.
cols());
265 const int M = A_rows;
266 const int N = B_cols;
267 const int K = A_cols;
269 assert(A_cols == B_rows &&
"op(A) cols must equal op(B) rows");
270 assert(
static_cast<int>(c->rows()) == M &&
static_cast<int>(c->cols()) == N &&
271 "Output matrix c must be sized MxN");
276 const int lda =
static_cast<int>(a.
cols());
277 const int ldb =
static_cast<int>(b.
cols());
278 const int ldc =
static_cast<int>(c->cols());
280 if constexpr (std::is_same_v<T, double>) {
281 cblas_dgemm(CblasRowMajor, ta, tb, M, N, K, 1.0, a.
data(), lda, b.
data(),
282 ldb, 0.0, c->data(), ldc);
284 }
else if constexpr (std::is_same_v<T, float>) {
285 cblas_sgemm(CblasRowMajor, ta, tb, M, N, K, 1.0f, a.
data(), lda, b.
data(),
286 ldb, 0.0f, c->data(), ldc);
288 }
else if constexpr (std::is_same_v<T, std::complex<double>>) {
289 const std::complex<double> alpha{1.0, 0.0};
290 const std::complex<double> beta{0.0, 0.0};
291 cblas_zgemm(CblasRowMajor, ta, tb, M, N, K, &alpha, a.
data(), lda, b.
data(),
292 ldb, &beta, c->data(), ldc);
294 }
else if constexpr (std::is_same_v<T, std::complex<float>>) {
295 const std::complex<float> alpha{1.0f, 0.0f};
296 const std::complex<float> beta{0.0f, 0.0f};
297 cblas_cgemm(CblasRowMajor, ta, tb, M, N, K, &alpha, a.
data(), lda, b.
data(),
298 ldb, &beta, c->data(), ldc);
301 static_assert(!
sizeof(T),
"GEMM: unsupported scalar type");
311 const auto N = A.rows();
312 assert(A.cols() == A.rows() &&
"Must be square");
324 for (std::size_t i = 0; i < N; ++i) {
325 for (std::size_t j = 0; j < N; ++j) {
326 const auto cij = C[i][j];
327 for (std::size_t b = 0; b < N; ++b) {
328 X[j][b] = cij * E[j][b];
332 for (std::size_t a = 0; a < N; ++a) {
333 for (std::size_t b = 0; b < N; ++b) {
334 M[a][b] += A[a][i] * Y[a][b] * D[i][b];
341template <
typename T,
bool PARALLEL>
345 const auto N = A.rows();
346 assert(A.cols() == A.rows() &&
"Must be square");
351 if constexpr (PARALLEL) {
353#pragma omp parallel for collapse(2)
354 for (std::size_t a = 0; a < N; ++a) {
355 for (std::size_t b = 0; b < N; ++b) {
356 const T *Ba = &B[a][0];
358 for (std::size_t i = 0; i < N; ++i) {
359 const auto AaiDib = A[a][i] * D[i][b];
360 const T *Ci = &C[i][0];
361 for (std::size_t j = 0; j < N; ++j) {
362 Mab += AaiDib * Ba[j] * Ci[j] * E[j][b];
371 for (std::size_t a = 0; a < N; ++a) {
372 const T *Ba = &B[a][0];
373 for (std::size_t b = 0; b < N; ++b) {
375 for (std::size_t i = 0; i < N; ++i) {
376 const auto AaiDib = A[a][i] * D[i][b];
377 const T *Ci = &C[i][0];
378 for (std::size_t j = 0; j < N; ++j) {
379 Mab += AaiDib * Ba[j] * Ci[j] * E[j][b];
391 if constexpr (std::is_same_v<T, double>) {
392 return gsl_matrix_view_array(m_data.data(), m_rows, m_cols);
393 }
else if constexpr (std::is_same_v<T, float>) {
394 return gsl_matrix_float_view_array(m_data.data(), m_rows, m_cols);
395 }
else if constexpr (std::is_same_v<T, std::complex<double>>) {
397 return gsl_matrix_complex_view_array(
398 reinterpret_cast<double *
>(m_data.data()), m_rows, m_cols);
399 }
else if constexpr (std::is_same_v<T, std::complex<float>>) {
400 return gsl_matrix_complex_float_view_array(
401 reinterpret_cast<float *
>(m_data.data()), m_rows, m_cols);
403 assert(
false &&
"as_gsl_view() only available for double/float (or complex "
410 if constexpr (std::is_same_v<T, double>) {
411 return gsl_matrix_const_view_array(m_data.data(), m_rows, m_cols);
412 }
else if constexpr (std::is_same_v<T, float>) {
413 return gsl_matrix_float_const_view_array(m_data.data(), m_rows, m_cols);
414 }
else if constexpr (std::is_same_v<T, std::complex<double>>) {
415 return gsl_matrix_complex_const_view_array(
416 reinterpret_cast<const double *
>(m_data.data()), m_rows, m_cols);
417 }
else if constexpr (std::is_same_v<T, std::complex<float>>) {
418 return gsl_matrix_complex_float_const_view_array(
419 reinterpret_cast<const float *
>(m_data.data()), m_rows, m_cols);
421 assert(
false &&
"as_gsl_view() only for available double/float (or complex "
428std::ostream &operator<<(std::ostream &os,
const Matrix<T> &a) {
429 for (
auto i = 0ul; i < a.rows(); ++i) {
430 for (
auto j = 0ul; j < a.cols(); ++j) {
431 os << a(i, j) <<
" ";
447 if constexpr (std::is_same_v<T, float> ||
448 std::is_same_v<T, std::complex<float>>) {
450 }
else if constexpr (std::is_same_v<T, double> ||
451 std::is_same_v<T, std::complex<double>>) {
466 for (
auto i = 0ul; i < lhs.
rows(); ++i) {
467 for (
auto j = 0ul; j < lhs.
cols(); ++j) {
469 if (std::abs(lhs(i, j) - rhs(i, j)) >
470 std::abs(eps * (lhs(i, j) + rhs(i, j))))
Row-major dense matrix with arithmetic and linear algebra support.
Definition Matrix.hpp:209
Matrix< T > & make_identity()
Constructs a diagonal unit matrix (identity), in place; only for square.
Definition Matrix.ipp:97
Matrix< T > & invert_in_place()
Inverts the matrix in place via LU decomposition (GSL).
Definition Matrix.ipp:37
T * data()
Pointer to first element; for std::complex<T> this is complex<T>*, not T*.
Definition Matrix.hpp:292
std::size_t rows() const
Return rows [major index size].
Definition Matrix.hpp:282
Matrix< T > conj() const
Returns conjugate of matrix.
Definition Matrix.ipp:117
Matrix< T > & operator+=(const Matrix< T > &rhs)
In-place elementwise addition; dimensions must match.
Definition Matrix.ipp:172
auto complex() const
Converts a real to complex matrix (changes type; returns a complex matrix)
Definition Matrix.ipp:159
auto imag() const
Returns imag part of complex matrix (changes type; returns a real matrix)
Definition Matrix.ipp:148
auto as_gsl_view()
Returns a GSL matrix view for use with GSL functions (no copy).
Definition Matrix.ipp:390
Matrix< T > & conj_in_place()
Conjugates matrix, in place.
Definition Matrix.ipp:128
Matrix< T > & operator/=(const T x)
In-place scalar divide: M_ij /= x.
Definition Matrix.ipp:194
Matrix< T > transpose() const
Returns the transpose of the matrix.
Definition Matrix.ipp:65
Matrix< T > & zero()
Sets all elements to zero, in place.
Definition Matrix.ipp:108
Matrix< T > & operator*=(const T x)
In-place scalar multiply: M_ij *= x.
Definition Matrix.ipp:188
std::size_t cols() const
Return columns [minor index size].
Definition Matrix.hpp:284
T determinant() const
Returns the determinant via LU decomposition (GSL).
Definition Matrix.ipp:9
Matrix< T > & mult_elements_by(const Matrix< T > &a)
Elementwise multiply in place: M_ij *= a_ij.
Definition Matrix.ipp:224
auto real() const
Returns real part of complex matrix (changes type; returns a real matrix)
Definition Matrix.ipp:137
Matrix< T > & operator-=(const Matrix< T > &rhs)
In-place elementwise subtraction; dimensions must match.
Definition Matrix.ipp:180
Linear algebra: matrices, vectors, views, and solvers.
Definition Matrix.hpp:54
void PENTA(const Matrix< T > &A, const Matrix< T > &B, const Matrix< T > &C, const Matrix< T > &D, const Matrix< T > &E, Matrix< T > *M)
5-matrix contraction for N*N matrices: M_ab = A_ai B_aj C_ij D_ib E_jb, without BLAS,...
Definition Matrix.ipp:342
void PENTA_GEMM(const Matrix< T > &A, const Matrix< T > &B, const Matrix< T > &C, const Matrix< T > &D, const Matrix< T > &E, Matrix< T > *M)
5-matrix contraction for N*N matrices: M_ab = A_ai B_aj C_ij D_ib E_jb, with BLAS
Definition Matrix.ipp:308
void GEMM(const Matrix< T > &a, const Matrix< T > &b, Matrix< T > *c, bool trans_A=false, bool trans_B=false)
Matrix multiplication C = op(A) * op(B) via CBLAS GEMM (row-major).
Definition Matrix.ipp:248
CBLAS_TRANSPOSE to_cblas_trans(bool trans)
Converts bool to CBLAS_TRANSPOSE enum (CblasTrans if true, CblasNoTrans if false)
Definition Matrix.hpp:547
constexpr auto myEps()
Default relative tolerance for equal(): 1e-6 for float, 1e-12 for double.
Definition Matrix.ipp:446
bool equal(const Matrix< T > &lhs, const Matrix< T > &rhs, T eps=myEps< T >())
Compares two matrices element-wise to within a relative tolerance.
Definition Matrix.ipp:461
Operator overloads for std::vector.
Definition Vector.hpp:503