19#ifndef THEORETICA_ALGEBRA_H
20#define THEORETICA_ALGEBRA_H
22#include "../complex/complex_types.h"
23#include "../core/core_traits.h"
24#include "../core/error.h"
43 template<
typename Matrix>
46 using Type = matrix_element_t<Matrix>;
48 for (
unsigned int i = 0; i < m.rows(); ++i)
49 for (
unsigned int j = 0; j < m.cols(); ++j)
50 m(i, j) = (Type) (i == j ?
nan() : 0);
61 template<
typename Vector>
64 using Type = vector_element_t<Vector>;
66 for (
unsigned int i = 0; i < v.size(); ++i)
78 template<
typename Vector>
92 template<
typename Matrix>
93 inline Matrix
make_error(
unsigned int rows,
unsigned int cols) {
95 result.resize(rows, cols);
103 template<
typename Matrix>
106 using Type = matrix_element_t<Matrix>;
108 for (
unsigned int i = 0; i < m.rows(); ++i)
109 for (
unsigned int j = 0; j < m.cols(); ++j)
110 m(i, j) = (Type) (i == j ? 1 : 0);
119 template<
typename Matrix>
122 using Type = matrix_element_t<Matrix>;
124 for (
unsigned int i = 0; i < m.rows(); ++i)
125 for (
unsigned int j = 0; j < m.cols(); ++j)
135 template<
typename Vector>
138 using Type = vector_element_t<Vector>;
140 for (
unsigned int i = 0; i < v.size(); ++i)
151 template<
typename Matrix1,
typename Matrix2>
152 inline Matrix1&
mat_copy(Matrix1& dest,
const Matrix2& src) {
154 dest.resize(src.rows(), src.cols());
156 for (
unsigned int i = 0; i < src.rows(); ++i)
157 for (
unsigned int j = 0; j < src.cols(); ++j)
158 dest(i, j) = src(i, j);
169 template<
typename Vector1,
typename Vector2>
170 inline Vector1&
vec_copy(Vector1& dest,
const Vector2& src) {
172 dest.resize(src.size());
174 for (
unsigned int i = 0; i < src.size(); ++i)
187 template<
typename Matrix>
188 inline Matrix&
mat_swap_rows(Matrix& A,
unsigned int row1,
unsigned int row2) {
190 using Type = matrix_element_t<Matrix>;
192 if (row1 >= A.rows()) {
197 if (row2 >= A.rows()) {
205 for (
unsigned int j = 0; j < A.cols(); ++j) {
207 const Type tmp = A(row1, j);
209 A(row1, j) = A(row2, j);
223 template<
typename Matrix>
224 inline Matrix&
mat_swap_cols(Matrix& A,
unsigned int col1,
unsigned int col2) {
226 using Type = matrix_element_t<Matrix>;
228 if (col1 >= A.cols()) {
233 if (col2 >= A.cols()) {
241 for (
unsigned int i = 0; i < A.cols(); ++i) {
243 const Type tmp = A(i, col1);
245 A(i, col1) = A(i, col2);
259 template<
typename Matrix,
typename Type = matrix_element_t<Matrix>>
262 const unsigned int count =
min(A.rows(), A.cols());
264 for (
unsigned int i = 0; i < count; ++i)
278 template<
typename Type>
290 template<
typename Type>
301 template<
typename Vector>
307 for (
unsigned int i = 1; i < v.size(); ++i)
318 template<
typename Vector>
319 inline auto norm(
const Vector& v) {
328 template<
typename Vector>
335 const auto m =
norm(v);
343 for (
unsigned int i = 0; i < r.size(); ++i)
353 template<
typename Vector>
356 const auto m =
norm(v);
364 for (
unsigned int i = 0; i < v.size(); ++i)
377 template<
typename Vector1,
typename Vector2>
378 inline auto dot(
const Vector1& v,
const Vector2& w) {
380 if(v.size() != w.size()) {
382 return vector_element_t<Vector1>(
nan());
388 for (
unsigned int i = 1; i < v.size(); ++i)
399 template<
typename Vector1,
typename Vector2>
400 inline Vector1
cross(
const Vector1& v1,
const Vector2& v2) {
415 v3[0] = v1[1] * v2[2] - v1[2] * v2[1];
416 v3[1] = v1[2] * v2[0] - v1[0] * v2[2];
417 v3[2] = v1[0] * v2[1] - v1[1] * v2[0];
427 template<
typename Matrix,
typename MatrixT = Matrix>
431 res.resize(m.cols(), m.rows());
433 for (
unsigned int i = 0; i < m.rows(); ++i)
434 for (
unsigned int j = 0; j < m.cols(); ++j)
444 template<
typename Matrix>
447 if(m.rows() != m.cols()) {
452 for (
unsigned int i = 0; i < m.rows(); ++i) {
453 for (
unsigned int j = 0; j < i; ++j) {
455 const auto buff = m(i, j);
471 template<
typename Matrix1,
typename Matrix2>
472 inline Matrix1&
transpose(Matrix1& dest,
const Matrix2& src) {
476 if(src.rows() != dest.cols()) {
481 if(src.cols() != dest.rows()) {
487 for (
unsigned int i = 0; i < src.rows(); ++i)
488 for (
unsigned int j = 0; j < src.cols(); ++j)
489 dest(j, i) = src(i, j);
499 template<
typename Matrix,
typename MatrixT = Matrix>
503 res.resize(m.cols(), m.rows());
505 for (
unsigned int i = 0; i < m.rows(); ++i)
506 for (
unsigned int j = 0; j < m.cols(); ++j)
518 template<
typename Matrix>
521 if(m.rows() != m.cols()) {
526 for (
unsigned int i = 0; i < m.rows(); ++i) {
527 for (
unsigned int j = 0; j < i; ++j) {
529 const auto buff = m(i, j);
546 template<
typename Matrix1,
typename Matrix2>
547 inline Matrix1&
hermitian(Matrix1& dest,
const Matrix2& src) {
551 if(src.rows() != dest.cols()) {
556 if(src.cols() != dest.rows()) {
562 for (
unsigned int i = 0; i < src.rows(); ++i)
563 for (
unsigned int j = 0; j < src.cols(); ++j)
573 template<
typename Matrix>
574 inline auto trace(
const Matrix& m) {
577 const size_t n =
min(m.rows(), m.cols());
579 for (
unsigned int i = 1; i < n; ++i)
591 template<
typename Matrix>
595 const size_t n =
min(m.rows(), m.cols());
597 for (
unsigned int i = 1; i < n; ++i)
608 template<
typename Field,
typename Matrix>
611 for (
unsigned int i = 0; i < m.rows(); ++i)
612 for (
unsigned int j = 0; j < m.cols(); ++j)
625 template<
typename Field,
typename Matrix1,
typename Matrix2>
626 inline Matrix1&
mat_scalmul(Matrix1& dest, Field a,
const Matrix2& src) {
628 if(src.rows() != dest.rows()) {
633 if(src.cols() != dest.cols()) {
638 for (
unsigned int i = 0; i < src.rows(); ++i)
639 for (
unsigned int j = 0; j < src.cols(); ++j)
640 dest(i, j) = a * src(i, j);
650 template<
typename Field,
typename Vector>
653 for (
unsigned int i = 0; i < v.size(); ++i)
666 template<
typename Field,
typename Vector1,
typename Vector2>
667 inline Vector1&
vec_scalmul(Vector1& dest, Field a,
const Vector2& src) {
669 if(src.size() != dest.size()) {
674 for (
unsigned int i = 0; i < src.size(); ++i)
675 dest[i] = a * src[i];
690 template<
typename Matrix,
typename Vector>
693 if(v.size() != A.cols()) {
699 res.resize(v.size());
702 for (
unsigned int i = 0; i < A.rows(); ++i)
703 for (
unsigned int j = 0; j < A.cols(); ++j)
704 res[i] += A(i, j) * v[j];
716 template<
typename Matrix,
typename Vector>
717 inline Vector
transform(
const Matrix& A,
const Vector& v) {
720 res.resize(v.size());
722 if(v.size() != A.cols()) {
729 for (
unsigned int i = 0; i < A.rows(); ++i)
730 for (
unsigned int j = 0; j < A.cols(); ++j)
731 res[i] += A(i, j) * v[j];
744 template<
typename Matrix,
typename Vector1,
typename Vector2>
745 inline Vector1&
transform(Vector1& res,
const Matrix& A,
const Vector2& v) {
747 if(v.size() != A.cols()) {
752 if(res.size() != v.size()) {
759 for (
unsigned int i = 0; i < A.rows(); ++i)
760 for (
unsigned int j = 0; j < A.cols(); ++j)
761 res[i] += A(i, j) * v[j];
775 template<
typename Matrix1,
typename Matrix2>
776 inline Matrix2&
mat_sum(Matrix1& A,
const Matrix2& B) {
778 if(A.rows() != B.rows()) {
783 if(A.cols() != B.cols()) {
788 for (
unsigned int i = 0; i < A.rows(); ++i)
789 for (
unsigned int j = 0; j < A.cols(); ++j)
790 A(i, j) = A(i, j) + B(i, j);
801 template<
typename Matrix1,
typename Matrix2,
typename Matrix3>
802 inline Matrix1&
mat_sum(Matrix1& res,
const Matrix2& A,
const Matrix3& B) {
804 if(A.rows() != B.rows()) {
809 if(A.cols() != B.cols()) {
814 if(res.rows() != A.rows()) {
819 if(res.cols() != A.cols()) {
824 for (
unsigned int i = 0; i < A.rows(); ++i)
825 for (
unsigned int j = 0; j < A.cols(); ++j)
826 res(i, j) = A(i, j) + B(i, j);
837 template<
typename Matrix1,
typename Matrix2>
838 inline Matrix2&
mat_diff(Matrix1& A,
const Matrix2& B) {
840 if(A.rows() != B.rows()) {
845 if(A.cols() != B.cols()) {
850 for (
unsigned int i = 0; i < A.rows(); ++i)
851 for (
unsigned int j = 0; j < A.cols(); ++j)
852 A(i, j) = A(i, j) - B(i, j);
863 template<
typename Matrix1,
typename Matrix2,
typename Matrix3>
864 inline Matrix1&
mat_diff(Matrix1& res,
const Matrix2& A,
const Matrix3& B) {
866 if(A.rows() != B.rows()) {
871 if(A.cols() != B.cols()) {
876 if(res.rows() != A.rows()) {
881 if(res.cols() != A.cols()) {
886 for (
unsigned int i = 0; i < A.rows(); ++i)
887 for (
unsigned int j = 0; j < A.cols(); ++j)
888 res(i, j) = A(i, j) - B(i, j);
900 template<
typename Field1,
typename Matrix1,
typename Field2,
typename Matrix2>
902 Field1 alpha, Matrix1& A, Field2 beta,
const Matrix2& B) {
904 if(A.rows() != B.rows()) {
909 if(A.cols() != B.cols()) {
914 for (
unsigned int i = 0; i < A.rows(); ++i)
915 for (
unsigned int j = 0; j < A.cols(); ++j)
916 A(i, j) = A(i, j) * alpha + B(i, j) * beta;
931 template<
typename Matrix1,
typename Field1,
typename Matrix2,
932 typename Field2,
typename Matrix3>
934 Matrix1& res, Field1 alpha,
const Matrix2& A, Field2 beta,
const Matrix3& B) {
936 if(A.rows() != B.rows()) {
941 if(A.cols() != B.cols()) {
946 if(res.rows() != A.rows()) {
952 if(res.cols() != A.cols()) {
957 for (
unsigned int i = 0; i < A.rows(); ++i)
958 for (
unsigned int j = 0; j < A.cols(); ++j)
959 res(i, j) = A(i, j) * alpha + B(i, j) * beta;
971 template<
typename Matrix1,
typename Matrix2,
typename Matrix3 = Matrix1>
972 inline Matrix3
mat_mul(
const Matrix1& A,
const Matrix2& B) {
975 R.resize(A.rows(), B.cols());
977 if(A.cols() != B.rows()) {
984 for (
unsigned int i = 0; i < A.rows(); ++i)
985 for (
unsigned int j = 0; j < B.cols(); ++j)
986 for (
unsigned int k = 0; k < A.cols(); ++k)
987 R(i, j) += A(i, k) * B(k, j);
1000 template<
typename Matrix1,
typename Matrix2,
typename Matrix3>
1001 inline Matrix1&
mat_mul(Matrix1& R,
const Matrix2& A,
const Matrix3& B) {
1003 if(R.rows() != A.rows()) {
1008 if(R.cols() != B.cols()) {
1013 if(A.cols() != B.rows()) {
1020 for (
unsigned int i = 0; i < A.rows(); ++i)
1021 for (
unsigned int j = 0; j < B.cols(); ++j)
1022 for (
unsigned int k = 0; k < A.cols(); ++k)
1023 R(i, j) += A(i, k) * B(k, j);
1039 template<
typename Matrix1,
typename Matrix2,
typename Matrix3 = Matrix1>
1043 R.resize(A.cols(), B.cols());
1045 if(A.rows() != B.rows()) {
1052 for (
unsigned int i = 0; i < R.rows(); ++i)
1053 for (
unsigned int j = 0; j < R.cols(); ++j)
1054 for (
unsigned int k = 0; k < A.rows(); ++k)
1055 R(i, j) += A(k, i) * B(k, j);
1071 template<
typename Matrix1,
typename Matrix2,
typename Matrix3 = Matrix1>
1075 R.resize(A.rows(), B.rows());
1077 if(A.cols() != B.cols()) {
1084 for (
unsigned int i = 0; i < R.rows(); ++i)
1085 for (
unsigned int j = 0; j < R.cols(); ++j)
1086 for (
unsigned int k = 0; k < A.cols(); ++k)
1087 R(i, j) += A(i, k) * B(j, k);
1099 template<
typename Matrix1,
typename Matrix2>
1103 if(A.rows() != B.rows() || A.cols() != B.cols())
1106 for (
unsigned int i = 0; i < A.rows(); ++i) {
1107 for (
unsigned int j = 0; j < A.cols(); ++j) {
1109 const auto diff =
abs(A(i, j) - B(i, j));
1111 if(diff > tolerance ||
is_nan(diff))
1125 template<
typename Vector1,
typename Vector2>
1126 inline Vector2&
vec_sum(Vector1& v1,
const Vector2& v2) {
1128 if(v1.size() != v2.size()) {
1133 for (
unsigned int i = 0; i < v1.size(); ++i)
1134 v1[i] = v1[i] + v2[i];
1145 template<
typename Vector1,
typename Vector2,
typename Vector3>
1146 inline Vector1&
vec_sum(Vector1& res,
const Vector2& v1,
const Vector3& v2) {
1148 if(v1.size() != v2.size()) {
1153 if(res.size() != v1.size()) {
1158 for (
unsigned int i = 0; i < v1.size(); ++i)
1159 res[i] = v1[i] + v2[i];
1170 template<
typename Vector1,
typename Vector2>
1171 inline Vector2&
vec_diff(Vector1& v1,
const Vector2& v2) {
1173 if(v1.size() != v2.size()) {
1178 for (
unsigned int i = 0; i < v1.size(); ++i)
1179 v1[i] = v1[i] - v2[i];
1190 template<
typename Vector1,
typename Vector2,
typename Vector3>
1191 inline Vector1&
vec_diff(Vector1& res,
const Vector2& v1,
const Vector3& v2) {
1193 if(v1.size() != v2.size()) {
1198 if(res.size() != v1.size()) {
1203 for (
unsigned int i = 0; i < v1.size(); ++i)
1204 res[i] = v1[i] - v2[i];
1216 template<
typename Matrix>
1218 return (m.rows() == m.cols());
1227 template<
typename Matrix>
1230 for (
unsigned int i = 0; i < m.rows(); ++i)
1231 for (
unsigned int j = 0; j < m.cols(); ++j)
1232 if(i != j &&
abs(m(i, j)) > tolerance)
1243 template<
typename Matrix>
1249 for (
unsigned int i = 0; i < m.rows(); ++i)
1250 for (
unsigned int j = 0; j < m.cols(); ++j)
1251 if(
abs(m(i, j) - m(j, i)) > tolerance)
1263 template<
typename Matrix>
1269 for (
unsigned int i = 0; i < m.rows(); ++i)
1270 for (
unsigned int j = i + 1; j < m.cols(); ++j)
1271 if (
abs(m(i, j)) > tolerance)
1283 template<
typename Matrix>
1289 for (
unsigned int i = 0; i < m.rows(); ++i)
1290 for (
unsigned int j = 0; j < i; ++j)
1291 if (
abs(m(i, j)) > tolerance)
1307 template<
typename Matrix, enable_matrix<Matrix> = true>
1310 long unsigned int n = 0;
1312 for (
unsigned int i = 0; i < A.rows(); ++i)
1313 for (
unsigned int j = 0; j < A.cols(); ++j)
1314 if (
abs(A(i, j)) > tolerance)
1317 return (
real(n) / A.rows()) / A.cols();
1330 template<
typename Matrix, enable_matrix<Matrix> = true>
1333 return 1.0 -
density(A, tolerance);
1347 template<
typename Matrix1,
typename Matrix2,
typename Matrix3>
1351 unsigned int err = 0;
1354 L.resize(A.rows(), A.cols());
1355 U.resize(A.rows(), A.cols());
1360 if (A.rows() != L.rows())
1363 if (A.cols() != L.cols())
1366 if (A.rows() != U.rows())
1369 if (A.cols() != U.cols())
1378 using Type = matrix_element_t<Matrix1>;
1381 for (
unsigned int i = 0; i < A.rows(); ++i)
1382 L(i, i) = (Type) 1.0;
1385 for(
unsigned int i = 0; i < A.rows(); ++i) {
1387 for(
unsigned int j = 0; j < A.rows(); ++j) {
1391 for(
unsigned int k = 0; k < i; ++k)
1395 for(
unsigned int j = i + 1; j < A.rows(); ++j) {
1399 for(
unsigned int k = 0; k < i; ++k)
1416 template<
typename Matrix>
1424 for(
unsigned int j = 0; j < A.cols(); ++j) {
1426 for(
unsigned int i = j + 1; i < A.rows(); ++i) {
1430 for(
unsigned int k = j + 1; k < A.rows(); ++k)
1445 template<
typename Matrix>
1449 L.resize(A.rows(), A.cols());
1450 using Type = matrix_element_t<Matrix>;
1464 for (
unsigned int i = 0; i < A.cols(); ++i) {
1466 for (
unsigned int j = 0; j <= i; ++j) {
1470 for (
unsigned int k = 1; k < j; ++k)
1475 const Type sqr_diag = A(j, j) -
sum;
1483 L(i, j) =
sqrt(sqr_diag);
1487 L(i, j) = (A(i, j) -
sum) / L(j, j);
1501 template<
typename Matrix>
1504 using Type = matrix_element_t<Matrix>;
1517 for (
unsigned int k = 0; k < A.rows(); ++k) {
1519 A(k, k) =
sqrt(A(k, k));
1521 for (
unsigned int i = k + 1; i < A.rows(); ++i)
1522 A(i, k) = A(i, k) / A(k, k);
1524 for (
unsigned int j = k + 1; j < A.rows(); ++j)
1525 for (
unsigned int i = j; i < A.cols(); ++i)
1530 for (
unsigned int i = 0; i < A.rows(); ++i)
1531 for (
unsigned int j = i + 1; j < A.rows(); ++j)
1532 A(i, j) = Type(0.0);
1547 template<
typename Matrix,
typename Vector>
1552 using Type = matrix_element_t<Matrix>;
1559 if (b.size() != L.rows()) {
1565 for (
unsigned int i = 0; i < L.cols(); ++i) {
1567 Type
sum = L(i, 0) * x[0];
1569 for (
unsigned int j = 1; j < i; ++j)
1570 sum += L(i, j) * x[j];
1577 x[i] = (b[i] -
sum) / L(i, i);
1589 template<
typename Matrix,
typename Vector>
1594 using Type = matrix_element_t<Matrix>;
1601 if (b.size() != U.rows()) {
1606 if (U.rows() == 1) {
1607 x[0] = b[0] / U(0, 0);
1612 for (
int i = U.rows() - 1; i >= 0; --i) {
1614 Type
sum = U(i, i + 1) * x[i + 1];
1616 for (
unsigned int j = i + 2; j < U.cols(); ++j)
1617 sum += U(i, j) * x[j];
1624 x[i] = (b[i] -
sum) / U(i, i);
1638 template<
typename Matrix,
typename Vector>
1648 err.resize(b.size());
1664 template<
typename Matrix,
typename Vector>
1672 if (A.rows() != b.size()) {
1677 using Type = matrix_element_t<Matrix>;
1680 for (
unsigned int i = 1; i < A.rows(); ++i) {
1682 Type
sum = Type(0.0);
1684 for (
unsigned int j = 0; j < i; ++j)
1685 sum += A(i, j) * b[j];
1691 for (
int i = A.rows() - 1; i >= 0; --i) {
1693 Type
sum = Type(0.0);
1695 for (
unsigned int j = i + 1; j < A.rows(); ++j)
1696 sum += A(i, j) * b[j];
1698 b[i] = (b[i] -
sum) / A(i, i);
1712 template<
typename Matrix,
typename Vector>
1733 template<
typename Matrix1,
typename Matrix2,
typename Vector>
1734 inline Vector
solve_lu(
const Matrix1& L,
const Matrix2& U,
const Vector& b) {
1748 if (L.rows() != U.rows()) {
1753 if (b.size() != L.rows()) {
1758 using Type1 = matrix_element_t<Matrix1>;
1761 for (
unsigned int i = 1; i < L.rows(); ++i) {
1763 Type1
sum = Type1(0.0);
1765 for (
unsigned int j = 0; j < i; ++j)
1766 sum += L(i, j) * x[j];
1771 using Type2 = matrix_element_t<Matrix2>;
1774 for (
int i = U.rows() - 1; i >= 0; --i) {
1776 Type2
sum = Type2(0.0);
1778 for (
unsigned int j = i + 1; j < U.rows(); ++j)
1779 sum += U(i, j) * x[j];
1786 x[i] = (x[i] -
sum) / U(i, i);
1804 template<
typename Matrix,
typename Vector>
1814 if (L.rows() != b.size()) {
1819 using Type = matrix_element_t<Matrix>;
1822 for (
unsigned int i = 0; i < L.rows(); ++i) {
1824 Type
sum = Type(0.0);
1826 for (
unsigned int j = 0; j < i; ++j)
1827 sum += L(i, j) * x[j];
1829 x[i] = (x[i] -
sum) / L(i, i);
1833 for (
int i = L.rows() - 1; i >= 0; --i) {
1835 Type
sum = Type(0.0);
1837 for (
unsigned int j = i + 1; j < L.rows(); ++j)
1838 sum += L(j, i) * x[j];
1840 x[i] = (x[i] -
sum) / L(i, i);
1853 template<
typename Matrix,
typename Vector>
1854 inline Vector
solve(
const Matrix& A,
const Vector& b) {
1869 template<
typename Matrix1,
typename Matrix2>
1870 inline Matrix1&
inverse(Matrix1& dest,
const Matrix2& src) {
1872 if(src.rows() != src.cols()) {
1877 if(dest.rows() != src.rows()) {
1882 if(dest.cols() != src.cols()) {
1887 using Type = matrix_element_t<Matrix2>;
1891 A.resize(src.rows(), src.cols());
1892 dest.resize(src.rows(), src.cols());
1897 for (
unsigned int i = 0; i < src.rows(); ++i) {
1901 if(A(i, i) == (Type) 0) {
1906 for (
unsigned int j = i + 1; j < src.rows(); ++j) {
1910 if(A(j, i) != (Type) 0) {
1912 for (
unsigned int k = 0; k < src.rows(); ++k) {
1914 dest(i, k) += dest(j, k);
1928 auto inv_pivot = ((Type) 1.0) / A(i, i);
1932 for (
unsigned int j = 0; j < src.rows(); ++j) {
1940 const auto coeff = A(j, i) * inv_pivot;
1942 for (
unsigned int k = 0; k < src.rows(); ++k) {
1943 A(j, k) -= coeff * A(i, k);
1944 dest(j, k) -= coeff * dest(i, k);
1949 for (
unsigned int j = 0; j < src.rows(); ++j) {
1950 A(i, j) *= inv_pivot;
1951 dest(i, j) *= inv_pivot;
1965 template<
typename Matrix,
typename MatrixInv = Matrix>
1968 res.resize(m.rows(), m.cols());
1979 template<
typename Matrix>
1982 if(m.rows() != m.cols()) {
1989 tmp.resize(m.rows(), m.cols());
2005 template<
typename Matrix>
2006 inline auto det(
const Matrix& A) {
2027 template<
typename Matrix,
typename Vector>
2030 const auto p =
dot(x, x);
2035 return vector_element_t<Vector>(
nan());
2053 template<
typename Matrix,
typename Vector>
2055 const Matrix& A,
const Vector& x,
2058 using Type = matrix_element_t<Matrix>;
2065 if (x.size() != A.rows()) {
2079 for (i = 1; i <= max_iter; ++i) {
2086 if (
norm(x_curr - x_prev) <= tolerance ||
2087 norm(x_curr + x_prev) <= tolerance)
2097 return dot(x_curr, A * x_curr);
2114 template<
typename Matrix,
typename Vector1,
typename Vector2 = Vector1>
2116 const Matrix& A,
const Vector1& x, Vector2& v,
2119 using Type = matrix_element_t<Matrix>;
2126 if (x.size() != A.rows()) {
2131 if (v.size() != x.size()) {
2145 for (i = 1; i <= max_iter; ++i) {
2152 if (
norm(x_curr - x_prev) <= tolerance ||
2153 norm(x_curr + x_prev) <= tolerance)
2166 return dot(x_curr, A * x_curr);
2181 template<
typename Matrix,
typename Vector>
2183 const Matrix& A,
const Vector& x,
2187 using Type = matrix_element_t<Matrix>;
2194 if (A.rows() != x.size()) {
2205 Vector x_curr = x_prev;
2211 for (i = 1; i <= max_iter; ++i) {
2219 if (
norm(x_curr - x_prev) <= tolerance ||
2220 norm(x_curr + x_prev) <= tolerance)
2231 return dot(x_curr, A * x_curr);
2248 template<
typename Matrix,
typename Vector1,
typename Vector2>
2250 const Matrix& A,
const Vector1& x, Vector2& v,
2253 using Type = matrix_element_t<Matrix>;
2260 if (A.rows() != x.size()) {
2271 Vector1 x_curr = x_prev;
2277 for (i = 1; i <= max_iter; ++i) {
2285 if (
norm(x_curr - x_prev) <= tolerance ||
2286 norm(x_curr + x_prev) <= tolerance)
2300 return dot(x_curr, A * x_curr);
2320 template<
typename Matrix,
typename Vector,
typename T = matrix_element_t<Matrix>>
2322 const Matrix& A,
const T& lambda,
const Vector& x,
2327 return make_error<Vector>(A.rows());
2330 if (A.rows() != x.size()) {
2332 return make_error<Vector>(A.rows());
2344 Vector v_curr = v_prev;
2350 for (i = 1; i <= max_iter; ++i) {
2358 if (
norm(v_curr - v_prev) <= tolerance ||
2359 norm(v_curr + v_prev) <= tolerance)
2386 template<
typename Matrix,
typename Vector,
typename T = matrix_element_t<Matrix>>
2388 const Matrix& A,
const T& lambda,
const Vector& x,
2391 using Type = matrix_element_t<Matrix>;
2398 if (A.rows() != x.size()) {
2407 Type lambda_prev = lambda;
2408 Type lambda_curr = lambda;
2415 for (i = 1; i <= max_iter; ++i) {
2419 lambda_prev = lambda_curr;
2433 if (
norm(x_curr - x_prev) <= tolerance ||
2434 norm(x_curr + x_prev) <= tolerance)
2462 template<
typename Matrix,
typename Vector1,
typename Vector2,
2463 typename T = matrix_element_t<Matrix>>
2465 const Matrix& A,
const T& lambda,
const Vector1& x, Vector2& v,
2468 using Type = matrix_element_t<Matrix>;
2475 if (A.rows() != x.size()) {
2484 Type lambda_prev = lambda;
2485 Type lambda_curr = lambda;
2492 for (i = 1; i <= max_iter; ++i) {
2496 lambda_prev = lambda_curr;
2510 if (
norm(x_curr - x_prev) <= tolerance ||
2511 norm(x_curr + x_prev) <= tolerance)
Complex number in algebraic form .
Definition complex.h:26
#define TH_MATH_ERROR(F_NAME, VALUE, EXCEPTION)
TH_MATH_ERROR is a macro which throws exceptions or modifies errno (depending on which compilation op...
Definition error.h:238
bool is_lower_triangular(const Matrix &m, real tolerance=ALGEBRA_ELEMENT_TOL)
Returns whether the matrix is lower triangular.
Definition algebra.h:1264
Vector1 cross(const Vector1 &v1, const Vector2 &v2)
Compute the cross product between two tridimensional vectors.
Definition algebra.h:400
Matrix2 & mat_sum(Matrix1 &A, const Matrix2 &B)
Sum two matrices and store the result in the first matrix.
Definition algebra.h:776
Matrix3 mat_transpose_mul(const Matrix1 &A, const Matrix2 &B)
Multiply the transpose of a matrix by another matrix, equivalent to the operation .
Definition algebra.h:1040
auto eigenpair_power(const Matrix &A, const Vector1 &x, Vector2 &v, real tolerance=ALGEBRA_EIGEN_TOL, unsigned int max_iter=ALGEBRA_EIGEN_ITER)
Find the biggest eigenvalue in module of a square matrix and its corresponding eigenvector (eigenpai...
Definition algebra.h:2115
Matrix & mat_scalmul(Field a, Matrix &m)
Multiply a matrix by a scalar of any compatible type.
Definition algebra.h:609
bool is_square(const Matrix &m)
Returns whether the matrix is square.
Definition algebra.h:1217
Matrix & invert(Matrix &m)
Invert the given matrix and overwrite it.
Definition algebra.h:1980
Vector make_error(unsigned int n)
Create a vector representing an error state, with all NaN values.
Definition algebra.h:79
auto eigenvalue_power(const Matrix &A, const Vector &x, real tolerance=ALGEBRA_EIGEN_TOL, unsigned int max_iter=ALGEBRA_EIGEN_ITER)
Find the biggest eigenvalue in module of a square matrix using the power method (Von Mises iteration...
Definition algebra.h:2054
auto trace(const Matrix &m)
Compute the trace of the given matrix.
Definition algebra.h:574
auto eigenpair_inverse(const Matrix &A, const Vector1 &x, Vector2 &v, real tolerance=ALGEBRA_EIGEN_TOL, unsigned int max_iter=ALGEBRA_EIGEN_ITER)
Find the eigenvalue with the smallest inverse of a square matrix and its corresponding eigenvector,...
Definition algebra.h:2249
Matrix & make_identity(Matrix &m)
Overwrite a matrix with the identity matrix.
Definition algebra.h:104
Matrix2 & mat_diff(Matrix1 &A, const Matrix2 &B)
Subtract two matrices and store the result in the first matrix.
Definition algebra.h:838
Matrix & mat_swap_cols(Matrix &A, unsigned int col1, unsigned int col2)
Swap two columns of a matrix, given the matrix and the two indices of the columns.
Definition algebra.h:224
Vector & make_normalized(Vector &v)
Normalize a given vector overwriting it.
Definition algebra.h:354
Vector1 & vec_copy(Vector1 &dest, const Vector2 &src)
Copy a vector by overwriting another.
Definition algebra.h:170
Vector & solve_lu_inplace(const Matrix &A, Vector &b)
Solve the linear system , finding , where the matrix A has already undergone in-place LU decompositio...
Definition algebra.h:1665
bool is_symmetric(const Matrix &m, real tolerance=ALGEBRA_ELEMENT_TOL)
Returns whether the matrix is symmetric.
Definition algebra.h:1244
Vector solve(const Matrix &A, const Vector &b)
Solve the linear system , finding using the best available algorithm.
Definition algebra.h:1854
Vector transform(const Matrix &A, const Vector &v)
Returns the matrix transformation of a vector.
Definition algebra.h:717
bool is_upper_triangular(const Matrix &m, real tolerance=ALGEBRA_ELEMENT_TOL)
Returns whether the matrix is upper triangular.
Definition algebra.h:1284
Vector & apply_transform(const Matrix &A, Vector &v)
Apply a matrix transformation to a vector and store the result in the vector.
Definition algebra.h:691
Vector eigenvector_inverse(const Matrix &A, const T &lambda, const Vector &x, real tolerance=ALGEBRA_EIGEN_TOL, unsigned int max_iter=ALGEBRA_EIGEN_ITER)
Find the eigenvector associated with a given approximated eigenvalue using the inverse power method.
Definition algebra.h:2321
void decompose_lu(const Matrix1 &A, Matrix2 &L, Matrix3 &U)
Decompose a square matrix to two triangular matrices, L and U where L is lower and U is upper,...
Definition algebra.h:1348
Matrix & mat_swap_rows(Matrix &A, unsigned int row1, unsigned int row2)
Swap two rows of a matrix, given the matrix and the two indices of the rows.
Definition algebra.h:188
Vector solve_lu(Matrix A, Vector b)
Solve the linear system , finding .
Definition algebra.h:1713
Vector solve_triangular(const Matrix &T, const Vector &b)
Solve the linear system for triangular .
Definition algebra.h:1639
Matrix & mat_error(Matrix &m)
Overwrite the given matrix with the error matrix with NaN values on the diagonal and zeroes everywher...
Definition algebra.h:44
Matrix decompose_cholesky(const Matrix &A)
Decompose a symmetric positive definite matrix into a triangular matrix so that using Cholesky decom...
Definition algebra.h:1446
auto dot(const Vector1 &v, const Vector2 &w)
Computes the dot product between two vectors, using the Hermitian form if needed.
Definition algebra.h:378
Matrix3 mat_mul(const Matrix1 &A, const Matrix2 &B)
Multiply two matrices and store the result in the first matrix, equivalent to the operation .
Definition algebra.h:972
Matrix & mat_zeroes(Matrix &m)
Overwrite a matrix with all zeroes.
Definition algebra.h:120
auto pair_inner_product(const Type &v_i, const Type &w_i)
Compute the contribution of the inner product between a pair of elements of two vectors,...
Definition algebra.h:279
auto diagonal_product(const Matrix &m)
Compute the product of the elements of the main diagonal of a generic matrix.
Definition algebra.h:592
auto eigenvalue_rayleigh(const Matrix &A, const T &lambda, const Vector &x, real tolerance=ALGEBRA_EIGEN_TOL, unsigned int max_iter=ALGEBRA_EIGEN_ITER)
Compute an eigenvalue of a square matrix using the Rayleigh quotient iteration method.
Definition algebra.h:2387
auto eigenpair_rayleigh(const Matrix &A, const T &lambda, const Vector1 &x, Vector2 &v, real tolerance=ALGEBRA_EIGEN_TOL, unsigned int max_iter=ALGEBRA_EIGEN_ITER)
Compute an eigenvalue of a square matrix and its corresponding eigenvector (eigenpair) using the Rayl...
Definition algebra.h:2464
auto rayleigh_quotient(const Matrix &A, const Vector &x)
Compute the Rayleigh quotient of a vector with respect to a matrix.
Definition algebra.h:2028
auto sqr_norm(const Vector &v)
Returns the square of the Euclidean/Hermitian norm of the given vector.
Definition algebra.h:302
auto det(const Matrix &A)
Compute the determinant of a square matrix.
Definition algebra.h:2006
Matrix & mat_shift_diagonal(Matrix &A, const Type &sigma)
Shift the diagonal of a matrix by the given amount, overwriting the matrix itself,...
Definition algebra.h:260
Vector solve_triangular_lower(const Matrix &L, const Vector &b)
Solve the linear system for lower triangular .
Definition algebra.h:1548
Matrix2 & mat_lincomb(Field1 alpha, Matrix1 &A, Field2 beta, const Matrix2 &B)
Compute the linear combination of two matrices and store the result in the first matrix.
Definition algebra.h:901
Matrix & decompose_cholesky_inplace(Matrix &A)
Decompose a symmetric positive definite matrix in-place, overwriting the starting matrix,...
Definition algebra.h:1502
Vector2 & vec_sum(Vector1 &v1, const Vector2 &v2)
Sum two vectors and store the result in the first vector.
Definition algebra.h:1126
Matrix3 mat_mul_transpose(const Matrix1 &A, const Matrix2 &B)
Multiply a matrix by the transpose of another matrix, equivalent to the operation .
Definition algebra.h:1072
Matrix1 & mat_copy(Matrix1 &dest, const Matrix2 &src)
Copy a matrix by overwriting another.
Definition algebra.h:152
Vector2 & vec_diff(Vector1 &v1, const Vector2 &v2)
Subtract two vectors and store the result in the first vector.
Definition algebra.h:1171
MatrixT hermitian(const Matrix &m)
Returns the hermitian of the given matrix.
Definition algebra.h:500
Matrix & make_transposed(Matrix &m)
Transpose the given matrix.
Definition algebra.h:445
Vector solve_triangular_upper(const Matrix &U, const Vector &b)
Solve the linear system for upper triangular .
Definition algebra.h:1590
Vector & vec_zeroes(Vector &v)
Overwrite a vector with all zeroes.
Definition algebra.h:136
Vector & vec_error(Vector &v)
Overwrite the given vector with the error vector with NaN values.
Definition algebra.h:62
real density(const Matrix &A, real tolerance=ALGEBRA_ELEMENT_TOL)
Compute the density of a matrix, counting the proportion of non-zero (bigger in module than the given...
Definition algebra.h:1308
Vector solve_cholesky(const Matrix &L, const Vector &b)
Solve a linear system defined by a symmetric positive definite matrix, using the Cholesky decomposit...
Definition algebra.h:1805
Matrix1 & inverse(Matrix1 &dest, const Matrix2 &src)
Invert the given matrix.
Definition algebra.h:1870
auto norm(const Vector &v)
Returns the Euclidean/Hermitian norm of the given vector.
Definition algebra.h:319
bool is_diagonal(const Matrix &m, real tolerance=ALGEBRA_ELEMENT_TOL)
Returns whether the matrix is diagonal.
Definition algebra.h:1228
bool mat_equals(const Matrix1 &A, const Matrix2 &B, real tolerance=ALGEBRA_ELEMENT_TOL)
Checks whether two matrices are equal.
Definition algebra.h:1100
Vector normalize(const Vector &v)
Returns the normalized vector.
Definition algebra.h:329
auto eigenvalue_inverse(const Matrix &A, const Vector &x, real tolerance=ALGEBRA_EIGEN_TOL, unsigned int max_iter=ALGEBRA_EIGEN_ITER)
Find the eigenvalue with the smallest inverse of a square matrix, using the inverse power method wit...
Definition algebra.h:2182
real sparsity(const Matrix &A, real tolerance=ALGEBRA_ELEMENT_TOL)
Compute the sparsity of a matrix, counting the proportion of zero (smaller in module than the given t...
Definition algebra.h:1331
Matrix & decompose_lu_inplace(Matrix &A)
Decompose a square matrix to two triangular matrices, L and U where L is lower and U is upper,...
Definition algebra.h:1417
Matrix & make_hermitian(Matrix &m)
Compute the hermitian of a given matrix and overwrite it.
Definition algebra.h:519
Vector & vec_scalmul(Field a, Vector &v)
Multiply a vector by a scalar of any compatible type.
Definition algebra.h:651
MatrixT transpose(const Matrix &m)
Returns the transpose of the given matrix.
Definition algebra.h:428
Main namespace of the library which contains all functions and objects.
Definition algebra.h:27
double real
A real number, defined as a floating point type.
Definition constants.h:207
auto min(const Vector &X)
Finds the minimum value inside a dataset.
Definition dataset.h:351
dual2 sqrt(dual2 x)
Compute the square root of a second order dual number.
Definition dual2_functions.h:54
bool is_nan(const T &x)
Check whether a generic variable is (equivalent to) a NaN number.
Definition error.h:94
dual2 abs(dual2 x)
Compute the absolute value of a second order dual number.
Definition dual2_functions.h:198
constexpr real ALGEBRA_ELEMENT_TOL
Tolerance for the elements of matrices.
Definition constants.h:279
auto sum(const Vector &X)
Compute the sum of a vector of real values using pairwise summation to reduce round-off error.
Definition dataset.h:219
dual2 conjugate(dual2 x)
Return the conjugate of a second order dual number.
Definition dual2_functions.h:35
TH_CONSTEXPR real nan()
Return a quiet NaN number in floating point representation.
Definition error.h:78
@ InvalidArgument
Invalid argument.
@ ImpossibleOperation
Mathematically impossible operation.
@ NoConvergence
Algorithm did not converge.
@ DivByZero
Division by zero.
constexpr real MACH_EPSILON
Machine epsilon for the real type.
Definition constants.h:216
constexpr real ALGEBRA_EIGEN_ITER
Maximum number of iterations for eigensolvers.
Definition constants.h:285
constexpr real ALGEBRA_EIGEN_TOL
Tolerance for eigensolvers.
Definition constants.h:282