19template <const
int n, const
int m,
typename T,
typename Mtype>
22template <const
int n, const
int m,
typename T =
double>
25template <
int n,
int m,
typename T>
28template <
int n,
typename T>
34template <
int n,
typename T>
40template <
int n,
typename T>
46template <
int n,
typename T>
62 static_assert(M::is_matrix() && M::rows() == M::columns(),
"SVD only for square matrix");
64 DiagonalMatrix<M::size(), hila::arithmetic_type<M>> singularvalues;
74 static_assert(M::is_matrix() && M::rows() == M::columns(),
75 "Eigenvalues only for square matrix");
105template <const
int n, const
int m,
typename T,
typename Mtype>
114 "Matrix requires Complex or arithmetic type");
117 using base_type = hila::arithmetic_type<T>;
118 using argument_type = T;
122 static constexpr bool is_matrix() {
133 return (n == 1 || m == 1);
154 template <
typename S,
int nn = n,
int mm = m,
155 std::enable_if_t<(hila::is_assignable<T &, S>::value && nn == mm),
int> = 0>
156 explicit inline Matrix_t(
const S rhs) {
157 for (
int i = 0; i < n; i++)
158 for (
int j = 0; j < n; j++) {
167 template <
typename S,
typename MT,
168 std::enable_if_t<hila::is_assignable<T &, S>::value,
int> = 0>
170 for (
int i = 0; i < n * m; i++) {
176 inline Matrix_t(
const std::nullptr_t &z) {
177 for (
int i = 0; i < n * m; i++) {
184 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value,
int> = 0>
185 inline Matrix_t(std::initializer_list<S> rhs) {
186 assert(rhs.size() == n * m &&
187 "Matrix/Vector initializer list size must match variable size");
189 for (
auto it = rhs.begin(); it != rhs.end(); it++, i++) {
196 template <
typename Tm =
Mtype,
197 std::enable_if_t<!std::is_same<Tm, Matrix<n, m, T>>::value,
int> = 0>
198 inline operator Mtype &() {
199 return *
reinterpret_cast<Mtype *
>(
this);
203 template <
typename Tm =
Mtype,
204 std::enable_if_t<!std::is_same<Tm, Matrix<n, m, T>>::value,
int> = 0>
205 inline operator const Mtype &()
const {
206 return *
reinterpret_cast<const Mtype *
>(
this);
247 template <
int q = n,
int p = m, std::enable_if_t<q == 1,
int> = 0>
252 template <
int q = n,
int p = m, std::enable_if_t<p == 1,
int> = 0>
253 static constexpr int size() {
257 template <
int q = n,
int p = m, std::enable_if_t<q == p,
int> = 0>
258 static constexpr int size() {
285#pragma hila loop_function
286 inline T
e(
const int i,
const int j)
const {
291 inline T &
e(
const int i,
const int j) const_function {
298#pragma hila loop_function
299 template <
int q = n,
int p = m, std::enable_if_t<(q == 1 || p == 1),
int> = 0>
300 inline T
e(
const int i)
const {
304 template <
int q = n,
int p = m, std::enable_if_t<(q == 1 || p == 1),
int> = 0>
305 inline T &
e(
const int i) const_function {
325#pragma hila loop_function
326 template <
int q = n,
int p = m, std::enable_if_t<(q == 1 || p == 1),
int> = 0>
331 template <
int q = n,
int p = m, std::enable_if_t<(q == 1 || p == 1),
int> = 0>
332 inline T &
operator[](
const int i) const_function {
344 for (
int i = 0; i < m; i++)
356 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value,
int> = 0>
358 for (
int i = 0; i < m; i++)
370 for (
int i = 0; i < n; i++)
394 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value,
int> = 0>
396 for (
int i = 0; i < n; i++)
407 static_assert(n == m,
"diagonal() method defined only for square matrices");
409 for (
int i = 0; i < n; i++)
410 res.
e(i) = (*this).e(i, i);
429 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value,
int> = 0>
431 static_assert(n == m,
"set_diagonal() method defined only for square matrices");
432 for (
int i = 0; i < n; i++)
433 (*this).e(i, i) = v.
e(i);
456 #pragma hila loop_function
457 template <
int mm = m, std::enable_if_t<mm == 1,
int> = 0>
462 template <
int mm = m, std::enable_if_t<mm == 1,
int> = 0>
495 for (
int i = 0; i < n * m; i++) {
521 template <
typename S,
int n2,
int m2>
523 if constexpr (n != n2 || m != m2)
526 for (
int i = 0; i < n; i++)
527 for (
int j = 0; j < m; j++) {
528 if (
e(i, j) != rhs.
e(i, j))
543 template <
typename S>
545 return !(*
this == rhs);
592 template <
typename S,
typename MT,
593 std::enable_if_t<hila::is_assignable<T &, S>::value,
int> = 0>
595 for (
int i = 0; i < n * m; i++) {
603 inline auto &
operator=(
const std::nullptr_t &z) out_only & {
604 for (
int i = 0; i < n * m; i++) {
612 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value && n == m,
int> = 0>
613 inline auto &
operator=(
const S rhs) out_only & {
615 for (
int i = 0; i < n; i++)
616 for (
int j = 0; j < n; j++) {
627 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value,
int> = 0>
629 static_assert(n == m,
630 "Assigning DiagonalMatrix to Matrix possible only for square matrices");
632 for (
int i = 0; i < n; i++)
633 for (
int j = 0; j < n; j++) {
645 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value,
int> = 0>
646 auto &
operator=(std::initializer_list<S> rhs) out_only & {
647 assert(rhs.size() == n * m &&
"Initializer list has a wrong size in assignment");
649 for (
auto it = rhs.begin(); it != rhs.end(); it++, i++) {
656 template <
typename S>
688 template <
typename S,
typename MT,
689 std::enable_if_t<hila::is_assignable<T &, S>::value,
int> = 0>
691 for (
int i = 0; i < n * m; i++) {
724 template <
typename S,
typename MT,
725 std::enable_if_t<hila::is_assignable<T &, S>::value,
int> = 0>
727 for (
int i = 0; i < n * m; i++) {
735 template <
typename S,
736 std::enable_if_t<hila::is_assignable<T &, hila::type_plus<T, S>>::value,
int> = 0>
739 static_assert(n == m,
"rows != columns : scalar addition possible for square matrix only!");
741 for (
int i = 0; i < n; i++) {
749 template <
typename S,
750 std::enable_if_t<hila::is_assignable<T &, hila::type_minus<T, S>>::value,
int> = 0>
752 static_assert(n == m,
753 "rows != columns : scalar subtraction possible for square matrix only!");
754 for (
int i = 0; i < n; i++) {
794 template <
int p,
typename S,
typename MT,
795 std::enable_if_t<hila::is_assignable<T &, hila::type_mul<T, S>>::value,
int> = 0>
797 static_assert(m == p,
"can't assign result of *= to lhs Matrix, because doing so "
798 "would change it's dimensions");
831 template <
typename S,
832 std::enable_if_t<hila::is_assignable<T &, hila::type_mul<T, S>>::value,
int> = 0>
834 for (
int i = 0; i < n * m; i++) {
861 template <
typename S,
862 std::enable_if_t<hila::is_assignable<T &, hila::type_div<T, S>>::value,
int> = 0>
864 for (
int i = 0; i < n * m; i++) {
875 template <
typename S,
876 std::enable_if_t<hila::is_assignable<T &, hila::type_plus<T, S>>::value,
int> = 0>
878 static_assert(n == m,
"Assigning DiagonalMatrix possible only for square matrix");
880 for (
int i = 0; i < n; i++)
885 template <
typename S,
886 std::enable_if_t<hila::is_assignable<T &, hila::type_plus<T, S>>::value,
int> = 0>
888 static_assert(n == m,
"Assigning DiagonalMatrix possible only for square matrix");
890 for (
int i = 0; i < n; i++)
898 template <
typename S,
899 std::enable_if_t<hila::is_assignable<T &, hila::type_mul<T, S>>::value,
int> = 0>
902 for (
int i = 0; i < n; i++)
903 for (
int j = 0; j < m; j++)
909 template <
typename S,
910 std::enable_if_t<hila::is_assignable<T &, hila::type_div<T, S>>::value,
int> = 0>
913 for (
int i = 0; i < n; i++)
914 for (
int j = 0; j < m; j++)
936 template <typename S, std::enable_if_t<hila::is_assignable<T &, S>::value,
int> = 0>
937 const auto &
fill(
const S rhs) out_only {
938 for (
int i = 0; i < n * m; i++)
951 template <
int mm = m,
952 typename Rtype =
typename std::conditional<n == m, Matrix_t, Matrix<m, n, T>>::type,
953 std::enable_if_t<(mm != 1 && n != 1),
int> = 0>
956 for (
int i = 0; i < n; i++)
957 for (
int j = 0; j < m; j++) {
958 res.e(j, i) =
e(i, j);
970 template <
int mm = m, std::enable_if_t<mm == 1,
int> = 0>
975 template <
int nn = n, std::enable_if_t<nn == 1 && m != 1,
int> = 0>
988 for (
int i = 0; i < n * m; i++) {
1001 template <
typename Rtype =
typename std::conditional<n == m, Mtype, Matrix<m, n, T>>::type>
1004 for (
int i = 0; i < n; i++)
1005 for (
int j = 0; j < m; j++) {
1006 res.e(j, i) =
::conj(
e(i, j));
1018 template <
typename Rtype =
typename std::conditional<n == m, Mtype, Matrix<m, n, T>>::type>
1033 for (
int i = 0; i < n * m; i++) {
1053 for (
int i = 0; i < m * n; i++) {
1066 for (
int i = 0; i < m * n; i++) {
1080 static_assert(n == m,
"trace not defined for non square matrices");
1082 for (
int i = 0; i < n; i++) {
1099 template <
int p,
int q,
typename S,
typename MT>
1102 static_assert(p == m && q == n,
"mul_trace(): argument matrix size mismatch");
1104 hila::type_mul<T, S> res = 0;
1105 for (
int i = 0; i < n; i++)
1106 for (
int j = 0; j < m; j++) {
1107 res +=
e(i, j) * rm.
e(j, i);
1118 hila::arithmetic_type<T> result(0);
1119 for (
int i = 0; i < n * m; i++) {
1131 template <
typename S = T,
1132 std::enable_if_t<hila::is_floating_point<hila::arithmetic_type<S>>::value,
int> = 0>
1133 hila::arithmetic_type<T>
norm()
const {
1137 template <
typename S = T,
1138 std::enable_if_t<!hila::is_floating_point<hila::arithmetic_type<S>>::value,
int> = 0>
1139 double norm()
const {
1147 template <typename S = T, std::enable_if_t<hila::is_arithmetic<S>::value,
int> = 0>
1150 for (
int i = 1; i < n * m; i++) {
1157 template <typename S = T, std::enable_if_t<hila::is_arithmetic<S>::value,
int> = 0>
1160 for (
int i = 1; i < n * m; i++) {
1167 auto max_abs()
const {
1168 hila::arithmetic_type<T> tres, res = 0;
1169 for (
int i = 0; i < n * m; i++) {
1178 auto min_abs()
const {
1179 hila::arithmetic_type<T> tres, res =
::abs(
c[0]);
1180 for (
int i = 1; i < n * m; i++) {
1219 template <
int p,
int q,
typename S,
typename R = hila::type_mul<T, S>>
1221 static_assert(m == 1 && q == 1 && p == n,
1222 "dot() product only for vectors of the same length");
1225 for (
int i = 0; i < n; i++) {
1262 template <
int p,
int q,
typename S,
typename R = hila::type_mul<T, S>>
1264 static_assert(m == 1 && q == 1,
"outer_product() only for vectors");
1267 for (
int i = 0; i < n; i++) {
1268 for (
int j = 0; j < p; j++) {
1298 static_assert(hila::is_floating_point<hila::arithmetic_type<T>>::value,
1299 "Matrix/Vector random() requires non-integral type elements");
1301 for (
int i = 0; i < n * m; i++) {
1316 static_assert(hila::is_floating_point<hila::arithmetic_type<T>>::value,
1317 "Matrix/Vector gaussian_random() requires non-integral type elements");
1320 if constexpr (hila::is_complex<T>::value) {
1321 for (
int i = 0; i < n * m; i++) {
1322 c[i].gaussian_random(width);
1329 for (
int i = 0; i < n * m - 1; i += 2) {
1331 c[i + 1] = gr * width;
1333 if constexpr ((n * m) % 2 > 0) {
1349 for (
int i = 0; i < m; i++)
1350 res.
set_column(i, this->column(permutation[i]));
1363 for (
int i = 0; i < n; i++)
1364 res.
set_row(i, this->row(permutation[i]));
1380 "permute() only for vectors, use permute_rows() or permute_columns() for matrices");
1381 static_assert(N ==
Mtype::size(),
"Incorrect size of permutation vector");
1384 for (
int i = 0; i < N; i++) {
1385 res[i] = (*this)[permutation[i]];
1414#pragma hila novector
1418 static_assert(n == 1 || m == 1,
"Sorting possible only for vectors");
1420 "Sorting possible only for arithmetic vector elements");
1421 static_assert(N ==
Mtype::size(),
"Incorrect size of permutation vector");
1423 for (
int i = 0; i < N; i++)
1425 if (hila::sort::unsorted == order) {
1429 if (hila::sort::ascending == order) {
1430 for (
int i = 1; i < N; i++) {
1431 for (
int j = i; j > 0 &&
c[permutation[j - 1]] >
c[permutation[j]]; j--)
1432 hila::swap(permutation[j], permutation[j - 1]);
1435 for (
int i = 1; i < N; i++) {
1436 for (
int j = i; j > 0 &&
c[permutation[j - 1]] <
c[permutation[j]]; j--)
1437 hila::swap(permutation[j], permutation[j - 1]);
1441 return this->
permute(permutation);
1444#pragma hila novector
1445 Mtype sort(hila::sort order = hila::sort::ascending)
const {
1446 static_assert(n == 1 || m == 1,
"Sorting possible only for vectors");
1449 return sort(permutation, order);
1472 template <
typename R,
typename Mt>
1474 static_assert(hila::contains_complex<T>::value || !hila::contains_complex<R>::value,
1475 "Cannot multiply real matrix with complex 2x2 block matrix");
1478 for (
int i = 0; i < m; i++) {
1479 a.
e(0) = this->
e(p, i);
1480 a.
e(1) = this->
e(q, i);
1484 this->
e(p, i) = a.e(0);
1485 this->
e(q, i) = a.e(1);
1501 template <
typename R,
typename Mt>
1503 static_assert(hila::contains_complex<T>::value || !hila::contains_complex<R>::value,
1504 "Cannot multiply real matrix with complex 2x2 block matrix");
1507 for (
int i = 0; i < n; i++) {
1508 a.
e(0) = this->
e(i, p);
1509 a.
e(1) = this->
e(i, q);
1513 this->
e(i, p) = a.
e(0);
1514 this->
e(i, q) = a.
e(1);
1527#pragma hila novector
1529#pragma hila novector
1531#pragma hila novector
1544#pragma hila novector
1545 template <
typename Et,
typename Mt,
typename MT>
1548 enum hila::sort sorted = hila::sort::unsorted)
const;
1553 template <
typename Et,
typename Mt,
typename MT>
1556 enum hila::sort sorted = hila::sort::unsorted)
const {
1560 eigenvaluevec = d.asArray().asVector();
1565#pragma hila novector
1566 template <
typename Et,
typename Mt,
typename MT>
1569 enum hila::sort sorted = hila::sort::unsorted)
const;
1574#pragma hila novector
1575 template <
typename Et,
typename Mt,
typename MT>
1578 enum hila::sort sorted = hila::sort::unsorted)
const;
1587 std::string
str(
int prec = 8,
char separator =
' ')
const {
1588 std::stringstream text;
1589 text.precision(prec);
1590 for (
int i = 0; i < n; i++) {
1591 for (
int j = 0; j < m; j++) {
1593 if (i < n - 1 || j < m - 1)
1619template <
int n,
int m,
typename T>
1625 using argument_type = T;
1678 Matrix &operator=(
const Matrix &v) out_only & =
default;
1691template <
typename T1,
typename T2>
1692struct matrix_op_type_s {
1693 using type =
Matrix<T1::rows(), T1::columns(), hila::ntype_op<T1, T2>>;
1697template <
typename T>
1698struct matrix_op_type_s<T, T> {
1703template <
typename T1,
typename T2>
1704using mat_x_mat_type =
typename matrix_op_type_s<T1, T2>::type;
1714template <
typename Mt,
typename S,
typename Enable =
void>
1715struct matrix_scalar_op_s {
1717 Matrix<Mt::rows(), Mt::columns(),
1721template <
typename Mt,
typename S>
1722struct matrix_scalar_op_s<
1724 typename std::enable_if_t<std::is_convertible<hila::type_plus<hila::number_type<Mt>, S>,
1725 hila::number_type<Mt>>::value>> {
1727 using type =
typename std::conditional<
1728 hila::is_floating_point<hila::arithmetic_type<Mt>>::value, Mt,
1729 Matrix<Mt::rows(), Mt::columns(),
1733template <
typename Mt,
typename S>
1734using mat_scalar_type =
typename matrix_scalar_op_s<Mt, S>::type;
1756template <
typename Mtype, std::enable_if_t<Mtype::is_matrix(),
int> = 0>
1758 return arg.transpose();
1776template <
typename Mtype, std::enable_if_t<Mtype::is_matrix(),
int> = 0>
1796template <
typename Mtype, std::enable_if_t<Mtype::is_matrix(),
int> = 0>
1798 return arg.adjoint();
1812template <
typename Mtype, std::enable_if_t<Mtype::is_matrix(),
int> = 0>
1814 return arg.adjoint();
1832template <
typename Mtype, std::enable_if_t<Mtype::is_matrix(),
int> = 0>
1853template <
typename Mtype, std::enable_if_t<Mtype::is_matrix(),
int> = 0>
1873template <
typename Mtype, std::enable_if_t<Mtype::is_matrix(),
int> = 0>
1893template <
typename Mtype, std::enable_if_t<Mtype::is_matrix(),
int> = 0>
1942template <
typename Mtype1,
typename Mtype2,
1943 std::enable_if_t<Mtype1::is_matrix() && Mtype2::is_matrix(),
int> = 0,
1944 typename Rtype = hila::mat_x_mat_type<Mtype1, Mtype2>,
1945 std::enable_if_t<!std::is_same<Mtype1, Rtype>::value,
int> = 0>
1948 constexpr int n = Mtype1::rows();
1949 constexpr int m = Mtype1::columns();
1951 static_assert(n == Mtype2::rows() && m == Mtype2::columns(),
"Matrix sizes do not match");
1954 for (
int i = 0; i < n; i++)
1955 for (
int j = 0; j < m; j++)
1956 r.e(i, j) = a.e(i, j) + b.e(i, j);
1962template <
typename Mtype1,
typename Mtype2,
1963 std::enable_if_t<Mtype1::is_matrix() && Mtype2::is_matrix(),
int> = 0,
1964 typename Rtype = hila::mat_x_mat_type<Mtype1, Mtype2>,
1965 std::enable_if_t<std::is_same<Mtype1, Rtype>::value,
int> = 0>
1966inline Rtype
operator+(Mtype1 a,
const Mtype2 &b) {
1968 constexpr int n = Mtype1::rows();
1969 constexpr int m = Mtype1::columns();
1971 static_assert(n == Mtype2::rows() && m == Mtype2::columns(),
"Matrix sizes do not match");
1973 for (
int i = 0; i < n; i++)
1974 for (
int j = 0; j < m; j++)
1975 a.e(i, j) += b.e(i, j);
2019template <
typename Mtype1,
typename Mtype2,
2020 std::enable_if_t<Mtype1::is_matrix() && Mtype2::is_matrix(),
int> = 0,
2021 typename Rtype = hila::mat_x_mat_type<Mtype1, Mtype2>>
2024 constexpr int n = Mtype1::rows();
2025 constexpr int m = Mtype1::columns();
2027 static_assert(n == Mtype2::rows() && m == Mtype2::columns(),
"Matrix sizes do not match");
2030 for (
int i = 0; i < n; i++)
2031 for (
int j = 0; j < m; j++)
2032 r.e(i, j) = a.e(i, j) - b.e(i, j);
2037template <
typename Mtype,
typename S,
2039 typename Rtype = hila::mat_scalar_type<Mtype, S>>
2043 "Matrix + scalar possible only for square matrix");
2046 for (
int i = 0; i < Rtype::rows(); i++)
2047 for (
int j = 0; j < Rtype::columns(); j++) {
2048 r.e(i, j) = a.
e(i, j);
2056template <
typename Mtype,
typename S,
2058 typename Rtype = hila::mat_scalar_type<Mtype, S>>
2062 "Matrix + scalar possible only for square matrix");
2067template <
typename Mtype,
typename S,
2069 typename Rtype = hila::mat_scalar_type<Mtype, S>>
2073 "Matrix - scalar possible only for square matrix");
2076 for (
int i = 0; i < Rtype::rows(); i++)
2077 for (
int j = 0; j < Rtype::columns(); j++) {
2078 r.e(i, j) = a.
e(i, j);
2086template <
typename Mtype,
typename S,
2088 typename Rtype = hila::mat_scalar_type<Mtype, S>>
2092 "Scalar - matrix possible only for square matrix");
2095 for (
int i = 0; i < Rtype::rows(); i++)
2096 for (
int j = 0; j < Rtype::columns(); j++) {
2097 r.e(i, j) = -a.
e(i, j);
2108template <
typename Mt, std::enable_if_t<Mt::is_matrix() && Mt::is_square(),
int> = 0>
2109inline Mt
operator*(
const Mt &a,
const Mt &b) {
2111 constexpr int n = Mt::rows();
2115 for (
int i = 0; i < n; i++)
2116 for (
int j = 0; j < n; j++) {
2117 res.e(i, j) = a.e(i, 0) * b.e(0, j);
2118 for (
int k = 1; k < n; k++) {
2119 res.e(i, j) += a.e(i, k) * b.e(k, j);
2182template <
typename Mt1,
typename Mt2,
2183 std::enable_if_t<Mt1::is_matrix() && Mt2::is_matrix() && !std::is_same<Mt1, Mt2>::value,
2185 typename R = hila::type_mul<hila::number_type<Mt1>, hila::number_type<Mt2>>,
2186 int n = Mt1::rows(),
int m = Mt2::columns()>
2189 constexpr int p = Mt1::columns();
2190 static_assert(p == Mt2::rows(),
"Matrix size: LHS columns != RHS rows");
2194 if constexpr (n > 1 && m > 1) {
2195 for (
int i = 0; i < n; i++)
2196 for (
int j = 0; j < m; j++) {
2197 res.
e(i, j) = a.e(i, 0) * b.e(0, j);
2198 for (
int k = 1; k < p; k++) {
2199 res.
e(i, j) += a.e(i, k) * b.e(k, j);
2202 }
else if constexpr (m == 1) {
2204 for (
int i = 0; i < n; i++) {
2205 res.
e(i) = a.e(i, 0) * b.e(0);
2206 for (
int k = 1; k < p; k++) {
2207 res.
e(i) += a.e(i, k) * b.e(k);
2210 }
else if constexpr (n == 1) {
2212 for (
int j = 0; j < m; j++) {
2213 res.
e(j) = a.e(0) * b.e(0, j);
2214 for (
int k = 1; k < p; k++) {
2215 res.
e(j) += a.e(k) * b.e(k, j);
2224template <
int m,
int n,
typename T1,
typename T2,
typename Rt = hila::ntype_op<T1, T2>>
2227 static_assert(m == n,
"Vector lengths do not match");
2230 res = A.
e(0) * B.
e(0);
2231 for (
int i = 1; i < m; i++) {
2232 res += A.
e(i) * B.
e(i);
2240template <
typename Mt,
typename S,
2242 typename Rt = hila::mat_scalar_type<Mt, S>>
2243inline Rt
operator*(
const Mt &mat,
const S &rhs) {
2245 for (
int i = 0; i < Rt::rows() * Rt::columns(); i++) {
2246 res.
c[i] = mat.c[i] * rhs;
2252template <
typename Mt,
typename S,
2254 typename Rt = hila::mat_scalar_type<Mt, S>>
2255inline Rt
operator*(
const S &rhs,
const Mt &mat) {
2282template <
typename Mt,
typename S,
2284 typename Rt = hila::mat_scalar_type<Mt, S>>
2285inline Rt
operator/(
const Mt &mat,
const S &rhs) {
2287 for (
int i = 0; i < Rt::rows() * Rt::columns(); i++) {
2288 res.c[i] = mat.c[i] / rhs;
2307template <
typename Mtype1,
typename Mtype2,
2308 std::enable_if_t<Mtype1::is_matrix() && Mtype2::is_matrix(),
int> = 0>
2311 static_assert(Mtype1::columns() == Mtype2::rows() && Mtype1::rows() == Mtype2::columns(),
2312 "mul_trace(a,b): matrix sizes are not compatible");
2313 return a.mul_trace(b);
2326template <
typename Mt, std::enable_if_t<Mt::is_matrix() && Mt::is_square(),
int> = 0>
2327inline void mult(
const Mt &a,
const Mt &b, out_only Mt &res) {
2328 constexpr int n = Mt::rows();
2330 for (i = 0; i < n; ++i) {
2331 for (j = 0; j < n; ++j) {
2332 res.e(i, j) = a.e(i, 0) * b.e(0, j);
2333 for (k = 1; k < n; ++k) {
2334 res.e(i, j) += a.e(i, k) * b.e(k, j);
2357template <
int n,
int m,
typename T,
typename MT>
2359 for (
int i = 0; i < n; i++)
2360 for (
int j = 0; j < m; j++) {
2362 if (i < n - 1 || j < m - 1)
2385template <
int n,
int m,
typename T,
typename MT>
2387 std::stringstream strm;
2388 strm.precision(prec);
2390 for (
int i = 0; i < n; i++)
2391 for (
int j = 0; j < m; j++) {
2393 if (i < n - 1 || j < m - 1)
2420template <
int n,
int m,
typename T,
typename MT>
2422 std::stringstream strm;
2423 strm.precision(prec);
2425 if constexpr (n == 1) {
2428 for (
int i = 0; i < n * m; i++)
2429 strm <<
' ' << prettyprint(A.
e(i), prec);
2436 std::vector<std::string> lines, columns;
2440 for (
int i = 0; i < n; i++)
2442 for (
int j = 0; j < m; j++) {
2444 for (
int i = 0; i < n; i++) {
2445 std::stringstream item;
2446 item << prettyprint(A.
e(i, j), prec);
2447 columns[i] = item.str();
2448 if (columns[i].size() > size)
2449 size = columns[i].size();
2451 for (
int i = 0; i < n; i++) {
2452 lines[i].append(size - columns[i].size(),
' ');
2453 lines[i].append(columns[i]);
2454 lines[i].append(1,
' ');
2457 for (
int i = 0; i < n - 1; i++) {
2458 strm << lines[i] <<
"]\n";
2460 strm << lines[n - 1] <<
"]";
2481template <
typename Mt, std::enable_if_t<Mt::is_matrix(),
int> = 0>
2483 return rhs.squarenorm();
2495template <
typename Mt, std::enable_if_t<Mt::is_matrix(),
int> = 0>
2504template <
typename Ntype,
typename T,
int n,
int m,
2505 std::enable_if_t<hila::is_arithmetic<T>::value,
int> = 0>
2508 for (
int i = 0; i < n * m; i++)
2509 res.
c[i] = mat.
c[i];
2513template <
typename Ntype,
typename T,
int n,
int m,
2514 std::enable_if_t<hila::is_complex<T>::value,
int> = 0>
2517 for (
int i = 0; i < n * m; i++)
2518 res.
c[i] = cast_to<Ntype>(mat.
c[i]);
2546template <
int n,
int m,
typename T,
typename MT>
2548 static_assert(n == m,
"exp() only for square matrices");
2550 hila::arithmetic_type<T> one = 1.0;
2554 for (
int k = order; k > 1; k--) {
2578template <
int n,
int m,
typename T,
typename MT>
2581 static_assert(n == m,
"chexp() only for square matrices");
2587 trpl[1] =
trace(pl[1]);
2589 for (i = 2; i <= n; ++i) {
2592 mult(pl[j], pl[j + k], pl[i]);
2593 trpl[i] =
trace(pl[i]);
2599 for (j = 1; j <= n; ++j) {
2601 for (i = 1; i <= j; ++i) {
2602 crpl[n - j] -= crpl[n - (j - i)] * trpl[i];
2610 hila::arithmetic_type<T>
2615 for (i = 0; i < n; ++i) {
2628 hila::arithmetic_type<T> s, rs = 1.0;
2630 for (j = n; j < mmax; ++j) {
2632 ch = -pal[n - 1] * crpl[0];
2637 for (i = 1; i < n; ++i) {
2638 ch = cho - pal[n - 1] * crpl[i];
2656 if (ttwpf == twpf) {
2667 for (i = 0; i < n; ++i) {
2668 for (j = 0; j < n; ++j) {
2670 omat.e(i, j) = al[0];
2674 for (k = 1; k < n; ++k) {
2675 omat.e(i, j) += al[k] * pl[k].e(i, j);
2682template <
int n,
int m,
typename T,
typename MT>
2685 static_assert(n == m,
"chexp() only for square matrices");
2686 chexp(mat, pl[0], pl);
2691template <
int n,
int m,
typename T,
typename MT>
2693 static_assert(n == m,
"chexp() only for square matrices");
2695 chexp(mat, omat, pl);
2700template <
int n,
int m,
typename T,
typename MT>
2702 static_assert(n == m,
"chexp() only for square matrices");
2704 chexp(mat, pl[0], pl);
2722template <
int n,
int m,
typename T,
typename MT>
2724 static_assert(n == m,
"chsexp() only for square matrices");
2737 for (k = 2; k <= n; ++k) {
2739 mult(mat, tB[1 - ip], tB[ip]);
2740 tc =
trace(tB[ip]) / k;
2748 hila::arithmetic_type<T> wpf = 1.0, twpf = 1.0, ttwpf;
2750 for (i = 0; i < n; ++i) {
2763 hila::arithmetic_type<T> s, rs = 1.0;
2765 for (j = n; j < mmax; ++j) {
2767 ch = -pal[n - 1] * crpl[0];
2772 for (i = 1; i < n; ++i) {
2773 ch = cho - pal[n - 1] * crpl[i];
2790 if (ttwpf == twpf) {
2803 tB[ip] *= al[n - 1];
2804 tB[ip] += al[n - 2];
2805 for (i = 2; i < n; ++i) {
2806 mult(tB[ip], mat, tB[1 - ip]);
2807 tB[1 - ip] += al[n - i - 1];
2817#include "datatypes/diagonal_matrix.h"
2819#include "datatypes/matrix_linalg.h"
Definition of Array class.
Array< n, m, T > sqrt(Array< n, m, T > a)
Square root.
Array< n, m, T > operator*(Array< n, m, T > a, const Array< n, m, T > &b)
Multiplication operator.
Array< n, m, Ntype > cast_to(const Array< n, m, T > &mat)
Array casting operation.
Define type DiagonalMatrix<n,T>
T e(const int i) const
Element access - e(i) gives diagonal element i.
The main matrix type template Matrix_t. This is a base class type for "useful" types which are deriv...
void set_column(int c, const Vector< n, S > &v)
get column of a matrix
void mult_by_2x2_left(int p, int q, const Matrix_t< 2, 2, R, Mt > &M)
Multiply -matrix from the left by matrix defined by sub matrix.
T trace() const
Computes Trace for Matrix.
Vector< n, T > column(int c) const
Returns column vector as value at index c.
hila::arithmetic_type< T > norm() const
Calculate vector norm - sqrt of squarenorm.
void set_row(int r, const RowVector< m, S > &v)
Set row of Matrix with RowVector if types are assignable.
auto & operator=(const Matrix_t< n, m, S, MT > &rhs) &
Assignment operator = to assign values to matrix.
static constexpr int rows()
Define constant methods rows(), columns() - may be useful in template code.
const auto & fill(const S rhs)
Matrix fill.
static constexpr int columns()
Returns column length.
Mtype & random()
dot with matrix - matrix
void mult_by_2x2_right(int p, int q, const Matrix_t< 2, 2, R, Mt > &M)
Multiply -matrix from the right by matrix defined by sub matrix.
Mtype sort(Vector< N, int > &permutation, hila::sort order=hila::sort::ascending) const
Sort method for Vector.
Matrix< n, m, hila::arithmetic_type< T > > imag() const
Returns imaginary part of Matrix or Vector.
void set_diagonal(const Vector< n, S > &v)
Set diagonal of square matrix to Vector which is passed to the method.
auto & operator-=(const Matrix_t< n, m, S, MT > &rhs) &
Subtract assign operator with matrix or scalar.
const RowVector< n, T > & transpose() const
Transpose of vector.
const auto & operator+() const
Unary + operator.
static constexpr int size()
Returns size of Vector or square Matrix.
Mtype permute(const Vector< N, int > &permutation) const
Permute elements of vector.
Matrix_t()=default
Define default constructors to ensure std::is_trivial.
Mtype conj() const
Returns complex conjugate of Matrix.
Mtype permute_rows(const Vector< n, int > &permutation) const
permute rows of Matrix
static constexpr bool is_vector()
Returns true if Matrix is a vector.
T det() const
determinant function - if matrix size is < 5, use Laplace, otherwise LU
auto abs() const
Returns absolute value of Matrix.
Matrix< n, p, R > outer_product(const Matrix< p, q, S > &rhs) const
Outer product.
T max() const
Find max or min value - only for arithmetic types.
Rtype dagger() const
Hermitian conjugate of matrix.
Matrix< n, m, hila::arithmetic_type< T > > real() const
Returns real part of Matrix or Vector.
bool operator!=(const Matrix< n, m, S > &rhs) const
Boolean operator != to check if matrices are exactly different.
int eigen_hermitean(DiagonalMatrix< n, Et > &eigenvalues, Matrix_t< n, n, Mt, MT > &eigenvectors, enum hila::sort sorted=hila::sort::unsorted) const
Calculate eigenvalues and -vectors of an hermitean (or symmetric) matrix.
int svd_pivot(Matrix_t< n, n, Mt, MT > &_U, DiagonalMatrix< n, Et > &_D, Matrix_t< n, n, Mt, MT > &_V, enum hila::sort sorted=hila::sort::unsorted) const
Singular value decomposition: divide matrix A = U S V*, where U,V unitary and S diagonal matrix of re...
bool operator==(const Matrix< n2, m2, S > &rhs) const
Boolean operator == to determine if two matrices are exactly the same.
T c[n *m]
The data as a one dimensional array.
auto & operator+=(const Matrix_t< n, m, S, MT > &rhs) &
Add assign operator with matrix or scalar.
Mtype & gaussian_random(base_type width=1.0)
Fills Matrix with gaussian random elements.
T e(const int i, const int j) const
Standard array indexing operation for matrices and vectors.
T det_laplace() const
Determinant of matrix, using Laplace method.
RowVector< m, T > row(int r) const
Return row of a matrix as a RowVector.
const DiagonalMatrix< n, T > & asDiagonalMatrix() const
Cast Vector to DiagonalMatrix.
Rtype transpose() const
Transpose of matrix.
auto & operator*=(const DiagonalMatrix< m, S > &rhs) &
mult and divide assign a diagonal - cols must match diagonal matrix rows
hila::type_mul< T, S > mul_trace(const Matrix_t< p, q, S, MT > &rm) const
Multiply with given matrix and compute trace of result.
T operator[](const int i) const
Indexing operation [] defined only for vectors.
int svd(Matrix_t< n, n, Mt, MT > &_U, DiagonalMatrix< n, Et > &_D, Matrix_t< n, n, Mt, MT > &_V, enum hila::sort sorted=hila::sort::unsorted) const
Singular value decomposition: divide matrix A = U S V*, where U,V unitary and S diagonal matrix of re...
auto & operator+=(const DiagonalMatrix< n, S > &rhs) &
add and sub assign a DiagonalMatrix
Mtype permute_columns(const Vector< m, int > &permutation) const
permute columns of Matrix
R dot(const Matrix< p, q, S > &rhs) const
Dot product.
T det_lu() const
following calculate the determinant of a square matrix det() is the generic interface,...
const Array< n, m, T > & asArray() const
Cast Matrix to Array.
Mtype operator-() const
Unary - operator.
auto & operator/=(const S rhs) &
Divide assign scalar.
DiagonalMatrix< n, T > diagonal()
Return diagonal of square matrix.
std::string str(int prec=8, char separator=' ') const
Rtype adjoint() const
Adjoint of matrix.
auto & operator*=(const Matrix_t< m, p, S, MT > &rhs) &
Multiply assign scalar or matrix.
hila::arithmetic_type< T > squarenorm() const
Calculate square norm - sum of squared elements.
static constexpr bool is_square()
Returns true if matrix is a square matrix.
Matrix class which defines matrix operations.
hila::arithmetic_type< T > base_type
std incantation for field types
Definition of Complex types.
T arg(const Complex< T > &a)
Return argument of Complex number.
This file defines all includes for HILA.
Matrix_t< n, m, T, MT > chsexp(const Matrix_t< n, m, T, MT > &mat)
Calculate exp of a square matrix.
void chexp(const Matrix_t< n, m, T, MT > &mat, Matrix_t< n, m, T, MT > &omat, Matrix_t< n, m, T, MT >(&pl)[n+1])
Calculate exp of a square matrix.
auto dagger(const Mtype &arg)
Return dagger of Matrix.
Rtype operator+(const Mtype1 &a, const Mtype2 &b)
Addition operator.
auto abs(const Mtype &arg)
Return absolute value Matrix or Vector.
Matrix_t< n, m, T, MT > exp(const Matrix_t< n, m, T, MT > &mat, const int order=20)
Calculate exp of a square matrix.
auto norm(const Mt &rhs)
Returns vector norm of Matrix.
auto conj(const Mtype &arg)
Return conjugate Matrix or Vector.
auto squarenorm(const Mt &rhs)
Returns square norm of Matrix.
auto transpose(const Mtype &arg)
Return transposed Matrix or Vector.
Rtype operator-(const Mtype1 &a, const Mtype2 &b)
Subtraction operator.
auto real(const Mtype &arg)
Return real of Matrix or Vector.
auto imag(const Mtype &arg)
Return imaginary of Matrix or Vector.
std::ostream & operator<<(std::ostream &strm, const Matrix_t< n, m, T, MT > &A)
Stream operator.
auto trace(const Mtype &arg)
Return trace of square Matrix.
auto mul_trace(const Mtype1 &a, const Mtype2 &b)
Returns trace of product of two matrices.
void mult(const Mt &a, const Mt &b, Mt &res)
compute product of two square matrices and write result to existing matrix
auto adjoint(const Mtype &arg)
Return adjoint Matrix.
Implement hila::swap for gauge fields.
double random()
Real valued uniform random number generator.
double gaussrand()
Gaussian random generation routine.
double gaussrand2(double &out2)
hila::gaussrand2 returns 2 Gaussian distributed random numbers with variance .
decltype(std::declval< A >()+std::declval< B >()) type_plus
std::string to_string(const Array< n, m, T > &A, int prec=8, char separator=' ')
Converts Array object to string.
type to store the return value of eigen_hermitean(): {E, U} where E is nxn DiagonalMatrix containing ...
hila::is_complex_or_arithmetic<T>::value
type to store the return combo of svd: {U, D, V} where U and V are nxn unitary / orthogonal,...