Theoretica
Scientific Computing
Loading...
Searching...
No Matches
algebra.h
Go to the documentation of this file.
1
18
19#ifndef THEORETICA_ALGEBRA_H
20#define THEORETICA_ALGEBRA_H
21
22#include "../complex/complex_types.h"
23#include "../core/core_traits.h"
24#include "../core/error.h"
25
26
27namespace theoretica {
28
29 namespace algebra {
30
31 // Error states for linear algebra types
32
39 template<typename Matrix>
40 inline Matrix& mat_error(Matrix& m) {
41
42 using Type = matrix_element_t<Matrix>;
43
44 for (unsigned int i = 0; i < m.rows(); ++i)
45 for (unsigned int j = 0; j < m.cols(); ++j)
46 m(i, j) = make_error<Type>();
47
48 return m;
49 }
50
51
57 template<typename Vector>
58 inline Vector& vec_error(Vector& v) {
59
60 using Type = vector_element_t<Vector>;
61
62 for (unsigned int i = 0; i < v.size(); ++i)
63 v[i] = make_error<Type>();
64
65 return v;
66 }
67 }
68
69
70 // These functions are outside the algebra namespace
71 // for compatibility with the generic make_error() function,
72 // defined in real_analysis.h
73
74
81 template <
82 typename Vector,
84 >
85 inline Vector make_error(unsigned int n) {
87 result.resize(n);
89 }
90
91
99 template <
100 typename Vector,
102 >
104
106 if (result.size() == 0)
107 result.resize(1);
108
110 }
111
112
120 template <
121 typename Matrix,
123 >
124 inline Matrix make_error(unsigned int rows, unsigned int cols) {
126 result.resize(rows, cols);
128 }
129
130
139 template <
140 typename Matrix,
142 >
144
146
147 if (result.rows() == 0 || result.cols() == 0)
148 result.resize(1, 1);
149
151 }
152
153
155 namespace algebra {
156
157
158 // Operations involving one matrix or vector
159
160
164 template<typename Matrix>
166
167 using Type = matrix_element_t<Matrix>;
168
169 for (unsigned int i = 0; i < m.rows(); ++i)
170 for (unsigned int j = 0; j < m.cols(); ++j)
171 m(i, j) = (Type) (i == j ? 1 : 0);
172
173 return m;
174 }
175
176
180 template<typename Matrix>
182
183 using Type = matrix_element_t<Matrix>;
184
185 for (unsigned int i = 0; i < m.rows(); ++i)
186 for (unsigned int j = 0; j < m.cols(); ++j)
187 m(i, j) = (Type) 0;
188
189 return m;
190 }
191
192
196 template<typename Vector>
198
199 using Type = vector_element_t<Vector>;
200
201 for (unsigned int i = 0; i < v.size(); ++i)
202 v[i] = Type(0.0);
203
204 return v;
205 }
206
207
212 template<typename Matrix1, typename Matrix2>
214
215 dest.resize(src.rows(), src.cols());
216
217 if (dest.rows() != src.rows()) {
218 TH_MATH_ERROR("algebra::mat_copy", dest.rows(), MathError::InvalidArgument);
219 return mat_error(dest);
220 }
221
222 if (dest.cols() != src.cols()) {
223 TH_MATH_ERROR("algebra::mat_copy", dest.cols(), MathError::InvalidArgument);
224 return mat_error(dest);
225 }
226
227 for (unsigned int i = 0; i < src.rows(); ++i)
228 for (unsigned int j = 0; j < src.cols(); ++j)
229 dest(i, j) = src(i, j);
230
231 return dest;
232 }
233
234
240 template<typename Vector1, typename Vector2>
242
243 dest.resize(src.size());
244
245 if (dest.size() != src.size()) {
246 TH_MATH_ERROR("algebra::vec_copy", dest.size(), MathError::InvalidArgument);
247 return vec_error(dest);
248 }
249
250 for (unsigned int i = 0; i < dest.size(); ++i)
251 dest[i] = src[i];
252
253 return dest;
254 }
255
256
263 template<typename Matrix>
264 inline Matrix& mat_swap_rows(Matrix& A, unsigned int row1, unsigned int row2) {
265
266 using Type = matrix_element_t<Matrix>;
267
268 if (row1 >= A.rows()) {
269 TH_MATH_ERROR("algebra::mat_swap_rows", row1, MathError::InvalidArgument);
270 return mat_error(A);
271 }
272
273 if (row2 >= A.rows()) {
274 TH_MATH_ERROR("algebra::mat_swap_rows", row2, MathError::InvalidArgument);
275 return mat_error(A);
276 }
277
278 if (row1 == row2)
279 return A;
280
281 for (unsigned int j = 0; j < A.cols(); ++j) {
282
283 const Type tmp = A(row1, j);
284
285 A(row1, j) = A(row2, j);
286 A(row2, j) = tmp;
287 }
288
289 return A;
290 }
291
292
299 template<typename Matrix>
300 inline Matrix& mat_swap_cols(Matrix& A, unsigned int col1, unsigned int col2) {
301
302 using Type = matrix_element_t<Matrix>;
303
304 if (col1 >= A.cols()) {
305 TH_MATH_ERROR("algebra::mat_swap_cols", col1, MathError::InvalidArgument);
306 return mat_error(A);
307 }
308
309 if (col2 >= A.cols()) {
310 TH_MATH_ERROR("algebra::mat_swap_cols", col2, MathError::InvalidArgument);
311 return mat_error(A);
312 }
313
314 if (col1 == col2)
315 return A;
316
317 for (unsigned int i = 0; i < A.cols(); ++i) {
318
319 const Type tmp = A(i, col1);
320
321 A(i, col1) = A(i, col2);
322 A(i, col2) = tmp;
323 }
324
325 return A;
326 }
327
328
335 template<typename Matrix, typename Type = matrix_element_t<Matrix>>
337
338 const unsigned int count = min(A.rows(), A.cols());
339
340 for (unsigned int i = 0; i < count; ++i)
341 A(i, i) += sigma;
342
343 return A;
344 }
345
346
354 template<typename Type>
355 inline auto pair_inner_product(const Type& v_i, const Type& w_i) {
356 return v_i * w_i;
357 }
358
366 template<typename Type>
368 return v_i * conjugate(w_i);
369 }
370
371
377 template<typename Vector>
378 inline auto sqr_norm(const Vector& v) {
379
380 auto sum = pair_inner_product(v[0], v[0]);
381
382 // Use conjugation for complex numbers
383 for (unsigned int i = 1; i < v.size(); ++i)
384 sum += pair_inner_product(v[i], v[i]);
385
386 return sum;
387 }
388
389
394 template<typename Vector>
395 inline auto norm(const Vector& v) {
396
398 }
399
400
404 template<typename Vector>
405 inline Vector normalize(const Vector& v) {
406
407 Vector r;
408 r.resize(v.size());
409 vec_copy(r, v);
410
411 const auto m = norm(v);
412
413 if(abs(m) < MACH_EPSILON) {
414 TH_MATH_ERROR("algebra::normalize", m, MathError::DivByZero);
415 vec_error(r);
416 return r;
417 }
418
419 for (unsigned int i = 0; i < r.size(); ++i)
420 r[i] /= m;
421
422 return r;
423 }
424
425
429 template<typename Vector>
431
432 const auto m = norm(v);
433
434 if(abs(m) < MACH_EPSILON) {
435 TH_MATH_ERROR("algebra::make_normalized", m, MathError::DivByZero);
436 vec_error(v);
437 return v;
438 }
439
440 for (unsigned int i = 0; i < v.size(); ++i)
441 v[i] /= m;
442
443 return v;
444 }
445
446
453 template<typename Vector1, typename Vector2>
454 inline auto dot(const Vector1& v, const Vector2& w) {
455
456 if(v.size() != w.size()) {
457 TH_MATH_ERROR("algebra::dot", v.size(), MathError::InvalidArgument);
459 }
460
461 auto sum = pair_inner_product(v[0], w[0]);
462
463 // Use conjugation for complex numbers
464 for (unsigned int i = 1; i < v.size(); ++i)
465 sum += pair_inner_product(v[i], w[i]);
466
467 return sum;
468 }
469
470
475 template<typename Vector1, typename Vector2>
476 inline Vector1 cross(const Vector1& v1, const Vector2& v2) {
477
478 Vector1 v3;
479 v3.resize(3);
480
481 if(v1.size() != 3) {
482 TH_MATH_ERROR("algebra::cross", v1.size(), MathError::InvalidArgument);
483 return vec_error(v3);
484 }
485
486 if(v2.size() != 3) {
487 TH_MATH_ERROR("algebra::cross", v2.size(), MathError::InvalidArgument);
488 return vec_error(v3);
489 }
490
491 v3[0] = v1[1] * v2[2] - v1[2] * v2[1];
492 v3[1] = v1[2] * v2[0] - v1[0] * v2[2];
493 v3[2] = v1[0] * v2[1] - v1[1] * v2[0];
494
495 return v3;
496 }
497
498
503 template<typename Matrix, typename MatrixT = Matrix>
504 inline MatrixT transpose(const Matrix& m) {
505
506 MatrixT res;
507 res.resize(m.cols(), m.rows());
508
509 for (unsigned int i = 0; i < m.rows(); ++i)
510 for (unsigned int j = 0; j < m.cols(); ++j)
511 res(j, i) = m(i, j);
512
513 return res;
514 }
515
516
520 template<typename Matrix>
522
523 if(m.rows() != m.cols()) {
524 TH_MATH_ERROR("algebra::make_transposed", m.rows(), MathError::InvalidArgument);
525 return mat_error(m);
526 }
527
528 for (unsigned int i = 0; i < m.rows(); ++i) {
529 for (unsigned int j = 0; j < i; ++j) {
530
531 const auto buff = m(i, j);
532 m(i, j) = m(j, i);
533 m(j, i) = buff;
534 }
535 }
536
537 return m;
538 }
539
540
547 template<typename Matrix1, typename Matrix2>
549
550 // Check that the two matrices have the correct
551 // number of rows and columns
552 if(src.rows() != dest.cols()) {
553 TH_MATH_ERROR("algebra::transpose", src.rows(), MathError::InvalidArgument);
554 return mat_error(dest);
555 }
556
557 if(src.cols() != dest.rows()) {
558 TH_MATH_ERROR("algebra::transpose", src.cols(), MathError::InvalidArgument);
559 return mat_error(dest);
560 }
561
562 // Overwrite dest with the transpose of src
563 for (unsigned int i = 0; i < src.rows(); ++i)
564 for (unsigned int j = 0; j < src.cols(); ++j)
565 dest(j, i) = src(i, j);
566
567 return dest;
568 }
569
570
575 template<typename Matrix, typename MatrixT = Matrix>
576 inline MatrixT hermitian(const Matrix& m) {
577
578 MatrixT res;
579 res.resize(m.cols(), m.rows());
580
581 for (unsigned int i = 0; i < m.rows(); ++i)
582 for (unsigned int j = 0; j < m.cols(); ++j)
583 res(j, i) = conjugate(m(i, j));
584
585 return res;
586 }
587
588
594 template<typename Matrix>
596
597 if(m.rows() != m.cols()) {
598 TH_MATH_ERROR("algebra::hermitian", m.rows(), MathError::InvalidArgument);
599 return mat_error(m);
600 }
601
602 for (unsigned int i = 0; i < m.rows(); ++i) {
603 for (unsigned int j = 0; j < i; ++j) {
604
605 const auto buff = m(i, j);
606 m(i, j) = conjugate(m(j, i));
607 m(j, i) = conjugate(buff);
608 }
609 }
610
611 return m;
612 }
613
614
622 template<typename Matrix1, typename Matrix2>
624
625 // Check that the two matrices have the correct
626 // number of rows and columns
627 if(src.rows() != dest.cols()) {
628 TH_MATH_ERROR("algebra::hermitian", src.rows(), MathError::InvalidArgument);
629 return mat_error(dest);
630 }
631
632 if(src.cols() != dest.rows()) {
633 TH_MATH_ERROR("algebra::hermitian", src.cols(), MathError::InvalidArgument);
634 return mat_error(dest);
635 }
636
637 // Overwrite dest with the transpose of src
638 for (unsigned int i = 0; i < src.rows(); ++i)
639 for (unsigned int j = 0; j < src.cols(); ++j)
640 dest(j, i) = conjugate(src(i, j));
641
642 return dest;
643 }
644
645
649 template<typename Matrix>
650 inline auto trace(const Matrix& m) {
651
652 auto sum = m(0, 0);
653 const size_t n = min(m.rows(), m.cols());
654
655 for (unsigned int i = 1; i < n; ++i)
656 sum += m(i, i);
657
658 return sum;
659 }
660
661
667 template<typename Matrix>
668 inline auto diagonal_product(const Matrix& m) {
669
670 auto mul = m(0, 0);
671 const size_t n = min(m.rows(), m.cols());
672
673 for (unsigned int i = 1; i < n; ++i)
674 mul *= m(i, i);
675
676 return mul;
677 }
678
679
684 template<typename Field, typename Matrix>
686
687 for (unsigned int i = 0; i < m.rows(); ++i)
688 for (unsigned int j = 0; j < m.cols(); ++j)
689 m(i, j) *= a;
690
691 return m;
692 }
693
694
701 template<typename Field, typename Matrix1, typename Matrix2>
703
704 if(src.rows() != dest.rows()) {
705 TH_MATH_ERROR("algebra::mat_scalmul", src.rows(), MathError::InvalidArgument);
706 return mat_error(dest);
707 }
708
709 if(src.cols() != dest.cols()) {
710 TH_MATH_ERROR("algebra::mat_scalmul", src.cols(), MathError::InvalidArgument);
711 return mat_error(dest);
712 }
713
714 for (unsigned int i = 0; i < src.rows(); ++i)
715 for (unsigned int j = 0; j < src.cols(); ++j)
716 dest(i, j) = a * src(i, j);
717
718 return dest;
719 }
720
721
726 template<typename Field, typename Vector>
728
729 for (unsigned int i = 0; i < v.size(); ++i)
730 v[i] *= a;
731
732 return v;
733 }
734
735
742 template<typename Field, typename Vector1, typename Vector2>
744
745 if(src.size() != dest.size()) {
746 TH_MATH_ERROR("algebra::vec_scalmul", src.size(), MathError::InvalidArgument);
747 return vec_error(dest);
748 }
749
750 for (unsigned int i = 0; i < src.size(); ++i)
751 dest[i] = a * src[i];
752
753 return dest;
754 }
755
756
757 // Operations involving a matrix and a vector
758
759
766 template<typename Matrix, typename Vector>
767 inline Vector& apply_transform(const Matrix& A, Vector& v) {
768
769 if(v.size() != A.cols()) {
770 TH_MATH_ERROR("algebra::apply_transform", v.size(), MathError::InvalidArgument);
771 return vec_error(v);
772 }
773
774 Vector res;
775 res.resize(A.rows());
777
778 for (unsigned int i = 0; i < A.rows(); ++i)
779 for (unsigned int j = 0; j < A.cols(); ++j)
780 res[i] += A(i, j) * v[j];
781
782 vec_copy(res, v);
783 return v;
784 }
785
786
798 template<typename Matrix, typename Vector>
799 inline Vector transform(const Matrix& A, const Vector& v) {
800
801 Vector res;
802 res.resize(A.rows());
803
804 if(v.size() != A.cols()) {
805 TH_MATH_ERROR("algebra::transform", v.size(), MathError::InvalidArgument);
806 return vec_error(res);
807 }
808
810
811 for (unsigned int i = 0; i < A.rows(); ++i)
812 for (unsigned int j = 0; j < A.cols(); ++j)
813 res[i] += A(i, j) * v[j];
814
815 return res;
816 }
817
818
826 template<typename Matrix, typename Vector1, typename Vector2>
827 inline Vector1& transform(Vector1& res, const Matrix& A, const Vector2& v) {
828
829 if(v.size() != A.cols()) {
830 TH_MATH_ERROR("algebra::transform", v.size(), MathError::InvalidArgument);
831 return vec_error(res);
832 }
833
834 if(res.size() != A.rows()) {
835 TH_MATH_ERROR("algebra::transform", res.size(), MathError::InvalidArgument);
836 return vec_error(res);
837 }
838
840
841 for (unsigned int i = 0; i < A.rows(); ++i)
842 for (unsigned int j = 0; j < A.cols(); ++j)
843 res[i] += A(i, j) * v[j];
844
845 return res;
846 }
847
848
849 // Operations involving multiple matrices or vectors
850
851
857 template<typename Matrix1, typename Matrix2>
858 inline Matrix2& mat_sum(Matrix1& A, const Matrix2& B) {
859
860 if(A.rows() != B.rows()) {
861 TH_MATH_ERROR("algebra::mat_sum", A.rows(), MathError::InvalidArgument);
862 return mat_error(A);
863 }
864
865 if(A.cols() != B.cols()) {
866 TH_MATH_ERROR("algebra::mat_sum", A.cols(), MathError::InvalidArgument);
867 return mat_error(A);
868 }
869
870 for (unsigned int i = 0; i < A.rows(); ++i)
871 for (unsigned int j = 0; j < A.cols(); ++j)
872 A(i, j) = A(i, j) + B(i, j);
873
874 return A;
875 }
876
877
883 template<typename Matrix1, typename Matrix2, typename Matrix3>
884 inline Matrix1& mat_sum(Matrix1& res, const Matrix2& A, const Matrix3& B) {
885
886 if(A.rows() != B.rows()) {
887 TH_MATH_ERROR("algebra::mat_sum", A.rows(), MathError::InvalidArgument);
888 return mat_error(res);
889 }
890
891 if(A.cols() != B.cols()) {
892 TH_MATH_ERROR("algebra::mat_sum", A.cols(), MathError::InvalidArgument);
893 return mat_error(res);
894 }
895
896 if(res.rows() != A.rows()) {
897 TH_MATH_ERROR("algebra::mat_sum", res.rows(), MathError::InvalidArgument);
898 return mat_error(res);
899 }
900
901 if(res.cols() != A.cols()) {
902 TH_MATH_ERROR("algebra::mat_sum", res.cols(), MathError::InvalidArgument);
903 return mat_error(res);
904 }
905
906 for (unsigned int i = 0; i < A.rows(); ++i)
907 for (unsigned int j = 0; j < A.cols(); ++j)
908 res(i, j) = A(i, j) + B(i, j);
909
910 return res;
911 }
912
913
919 template<typename Matrix1, typename Matrix2>
920 inline Matrix2& mat_diff(Matrix1& A, const Matrix2& B) {
921
922 if(A.rows() != B.rows()) {
923 TH_MATH_ERROR("algebra::mat_diff", A.rows(), MathError::InvalidArgument);
924 return mat_error(A);
925 }
926
927 if(A.cols() != B.cols()) {
928 TH_MATH_ERROR("algebra::mat_diff", A.cols(), MathError::InvalidArgument);
929 return mat_error(A);
930 }
931
932 for (unsigned int i = 0; i < A.rows(); ++i)
933 for (unsigned int j = 0; j < A.cols(); ++j)
934 A(i, j) = A(i, j) - B(i, j);
935
936 return A;
937 }
938
939
945 template<typename Matrix1, typename Matrix2, typename Matrix3>
946 inline Matrix1& mat_diff(Matrix1& res, const Matrix2& A, const Matrix3& B) {
947
948 if(A.rows() != B.rows()) {
949 TH_MATH_ERROR("algebra::mat_diff", A.rows(), MathError::InvalidArgument);
950 return mat_error(res);
951 }
952
953 if(A.cols() != B.cols()) {
954 TH_MATH_ERROR("algebra::mat_diff", A.cols(), MathError::InvalidArgument);
955 return mat_error(res);
956 }
957
958 if(res.rows() != A.rows()) {
959 TH_MATH_ERROR("algebra::mat_diff", res.rows(), MathError::InvalidArgument);
960 return mat_error(res);
961 }
962
963 if(res.cols() != A.cols()) {
964 TH_MATH_ERROR("algebra::mat_diff", res.cols(), MathError::InvalidArgument);
965 return mat_error(res);
966 }
967
968 for (unsigned int i = 0; i < A.rows(); ++i)
969 for (unsigned int j = 0; j < A.cols(); ++j)
970 res(i, j) = A(i, j) - B(i, j);
971
972 return res;
973 }
974
975
982 template<typename Field1, typename Matrix1, typename Field2, typename Matrix2>
984 Field1 alpha, Matrix1& A, Field2 beta, const Matrix2& B) {
985
986 if(A.rows() != B.rows()) {
987 TH_MATH_ERROR("algebra::mat_lincomb", A.rows(), MathError::InvalidArgument);
988 return mat_error(A);
989 }
990
991 if(A.cols() != B.cols()) {
992 TH_MATH_ERROR("algebra::mat_lincomb", A.cols(), MathError::InvalidArgument);
993 return mat_error(A);
994 }
995
996 for (unsigned int i = 0; i < A.rows(); ++i)
997 for (unsigned int j = 0; j < A.cols(); ++j)
998 A(i, j) = A(i, j) * alpha + B(i, j) * beta;
999
1000 return A;
1001 }
1002
1003
1013 template<typename Matrix1, typename Field1, typename Matrix2,
1014 typename Field2, typename Matrix3>
1016 Matrix1& res, Field1 alpha, const Matrix2& A, Field2 beta, const Matrix3& B) {
1017
1018 if(A.rows() != B.rows()) {
1019 TH_MATH_ERROR("algebra::mat_lincomb", A.rows(), MathError::InvalidArgument);
1020 return mat_error(res);
1021 }
1022
1023 if(A.cols() != B.cols()) {
1024 TH_MATH_ERROR("algebra::mat_lincomb", A.cols(), MathError::InvalidArgument);
1025 return mat_error(res);
1026 }
1027
1028 if(res.rows() != A.rows()) {
1029 TH_MATH_ERROR("algebra::mat_lincomb", res.rows(), MathError::InvalidArgument);
1030 return mat_error(res);
1031 return res;
1032 }
1033
1034 if(res.cols() != A.cols()) {
1035 TH_MATH_ERROR("algebra::mat_lincomb", res.cols(), MathError::InvalidArgument);
1036 return mat_error(res);
1037 }
1038
1039 for (unsigned int i = 0; i < A.rows(); ++i)
1040 for (unsigned int j = 0; j < A.cols(); ++j)
1041 res(i, j) = A(i, j) * alpha + B(i, j) * beta;
1042
1043 return A;
1044 }
1045
1046
1053 template<typename Matrix1, typename Matrix2, typename Matrix3 = Matrix1>
1054 inline Matrix3 mat_mul(const Matrix1& A, const Matrix2& B) {
1055
1056 Matrix3 R;
1057 R.resize(A.rows(), B.cols());
1058
1059 if(A.cols() != B.rows()) {
1060 TH_MATH_ERROR("algebra::mat_mul", A.cols(), MathError::InvalidArgument);
1061 return mat_error(R);
1062 }
1063
1064 mat_zeroes(R);
1065
1066 for (unsigned int i = 0; i < A.rows(); ++i)
1067 for (unsigned int j = 0; j < B.cols(); ++j)
1068 for (unsigned int k = 0; k < A.cols(); ++k)
1069 R(i, j) += A(i, k) * B(k, j);
1070
1071 return R;
1072 }
1073
1074
1082 template<typename Matrix1, typename Matrix2, typename Matrix3>
1083 inline Matrix1& mat_mul(Matrix1& R, const Matrix2& A, const Matrix3& B) {
1084
1085 if(R.rows() != A.rows()) {
1086 TH_MATH_ERROR("algebra::mat_mul", R.rows(), MathError::InvalidArgument);
1087 return mat_error(R);
1088 }
1089
1090 if(R.cols() != B.cols()) {
1091 TH_MATH_ERROR("algebra::mat_mul", R.cols(), MathError::InvalidArgument);
1092 return mat_error(R);
1093 }
1094
1095 if(A.cols() != B.rows()) {
1096 TH_MATH_ERROR("algebra::mat_mul", A.cols(), MathError::InvalidArgument);
1097 return mat_error(R);
1098 }
1099
1100 mat_zeroes(R);
1101
1102 for (unsigned int i = 0; i < A.rows(); ++i)
1103 for (unsigned int j = 0; j < B.cols(); ++j)
1104 for (unsigned int k = 0; k < A.cols(); ++k)
1105 R(i, j) += A(i, k) * B(k, j);
1106
1107 return R;
1108 }
1109
1110
1121 template<typename Matrix1, typename Matrix2, typename Matrix3 = Matrix1>
1122 inline Matrix3 mat_transpose_mul(const Matrix1& A, const Matrix2& B) {
1123
1124 Matrix3 R;
1125 R.resize(A.cols(), B.cols());
1126
1127 if(A.rows() != B.rows()) {
1128 TH_MATH_ERROR("algebra::mat_transpose_mul", A.rows(), MathError::InvalidArgument);
1129 return mat_error(R);
1130 }
1131
1132 mat_zeroes(R);
1133
1134 for (unsigned int i = 0; i < R.rows(); ++i)
1135 for (unsigned int j = 0; j < R.cols(); ++j)
1136 for (unsigned int k = 0; k < A.rows(); ++k)
1137 R(i, j) += A(k, i) * B(k, j);
1138
1139 return R;
1140 }
1141
1142
1153 template<typename Matrix1, typename Matrix2, typename Matrix3 = Matrix1>
1154 inline Matrix3 mat_mul_transpose(const Matrix1& A, const Matrix2& B) {
1155
1156 Matrix3 R;
1157 R.resize(A.rows(), B.rows());
1158
1159 if(A.cols() != B.cols()) {
1160 TH_MATH_ERROR("algebra::mat_mul_transpose", A.rows(), MathError::InvalidArgument);
1161 return mat_error(R);
1162 }
1163
1164 mat_zeroes(R);
1165
1166 for (unsigned int i = 0; i < R.rows(); ++i)
1167 for (unsigned int j = 0; j < R.cols(); ++j)
1168 for (unsigned int k = 0; k < A.cols(); ++k)
1169 R(i, j) += A(i, k) * B(j, k);
1170
1171 return R;
1172 }
1173
1174
1181 template<typename Matrix1, typename Matrix2>
1182 inline bool mat_equals(
1183 const Matrix1& A, const Matrix2& B, real tolerance = ALGEBRA_ELEMENT_TOL) {
1184
1185 if(A.rows() != B.rows() || A.cols() != B.cols())
1186 return false;
1187
1188 for (unsigned int i = 0; i < A.rows(); ++i) {
1189 for (unsigned int j = 0; j < A.cols(); ++j) {
1190
1191 const auto diff = abs(A(i, j) - B(i, j));
1192
1193 if(diff > tolerance || is_nan(diff))
1194 return false;
1195 }
1196 }
1197
1198 return true;
1199 }
1200
1201
1207 template<typename Vector1, typename Vector2>
1208 inline Vector2& vec_sum(Vector1& v1, const Vector2& v2) {
1209
1210 if(v1.size() != v2.size()) {
1211 TH_MATH_ERROR("algebra::vec_sum", v1.size(), MathError::InvalidArgument);
1212 return vec_error(v1);
1213 }
1214
1215 for (unsigned int i = 0; i < v1.size(); ++i)
1216 v1[i] = v1[i] + v2[i];
1217
1218 return v1;
1219 }
1220
1221
1227 template<typename Vector1, typename Vector2, typename Vector3>
1228 inline Vector1& vec_sum(Vector1& res, const Vector2& v1, const Vector3& v2) {
1229
1230 if(v1.size() != v2.size()) {
1231 TH_MATH_ERROR("algebra::vec_sum", v1.size(), MathError::InvalidArgument);
1232 return vec_error(res);
1233 }
1234
1235 if(res.size() != v1.size()) {
1236 TH_MATH_ERROR("algebra::vec_sum", res.size(), MathError::InvalidArgument);
1237 return vec_error(res);
1238 }
1239
1240 for (unsigned int i = 0; i < v1.size(); ++i)
1241 res[i] = v1[i] + v2[i];
1242
1243 return res;
1244 }
1245
1246
1252 template<typename Vector1, typename Vector2>
1253 inline Vector2& vec_diff(Vector1& v1, const Vector2& v2) {
1254
1255 if(v1.size() != v2.size()) {
1256 TH_MATH_ERROR("algebra::vec_diff", v1.size(), MathError::InvalidArgument);
1257 return vec_error(v1);
1258 }
1259
1260 for (unsigned int i = 0; i < v1.size(); ++i)
1261 v1[i] = v1[i] - v2[i];
1262
1263 return v1;
1264 }
1265
1266
1272 template<typename Vector1, typename Vector2, typename Vector3>
1273 inline Vector1& vec_diff(Vector1& res, const Vector2& v1, const Vector3& v2) {
1274
1275 if(v1.size() != v2.size()) {
1276 TH_MATH_ERROR("algebra::vec_diff", v1.size(), MathError::InvalidArgument);
1277 return vec_error(res);
1278 }
1279
1280 if(res.size() != v1.size()) {
1281 TH_MATH_ERROR("algebra::vec_diff", res.size(), MathError::InvalidArgument);
1282 return vec_error(res);
1283 }
1284
1285 for (unsigned int i = 0; i < v1.size(); ++i)
1286 res[i] = v1[i] - v2[i];
1287
1288 return res;
1289 }
1290
1291
1297 template<typename ReturnVector, typename Vector, typename Matrix>
1298 inline ReturnVector vec_mat_mul(const Vector& v, const Matrix& A) {
1299
1301 res.resize(A.cols());
1302
1303 if(v.size() != A.rows()) {
1304 TH_MATH_ERROR("algebra::vec_mul", v.size(), MathError::InvalidArgument);
1305 return vec_error(res);
1306 }
1307
1308 if(res.size() != A.cols()) {
1309 TH_MATH_ERROR("algebra::vec_mul", res.size(), MathError::ImpossibleOperation);
1310 return vec_error(res);
1311 }
1312
1313 vec_zeroes(res);
1314
1315 for (unsigned int i = 0; i < A.rows(); ++i)
1316 for (unsigned int j = 0; j < A.cols(); ++j)
1317 res[j] += v[i] * A(i, j);
1318
1319 return res;
1320 }
1321
1322
1329 template<typename Vector1, typename Vector2, typename Matrix>
1330 inline Vector1& vec_mat_mul(Vector1& res, const Vector2& v, const Matrix& A) {
1331
1332 if(v.size() != A.rows()) {
1333 TH_MATH_ERROR("algebra::vec_mul", v.size(), MathError::InvalidArgument);
1334 return vec_error(res);
1335 }
1336
1337 if(res.size() != A.cols()) {
1338 TH_MATH_ERROR("algebra::vec_mul", res.size(), MathError::InvalidArgument);
1339 return vec_error(res);
1340 }
1341
1342 vec_zeroes(res);
1343
1344 for (unsigned int i = 0; i < A.rows(); ++i)
1345 for (unsigned int j = 0; j < A.cols(); ++j)
1346 res[j] += v[i] * A(i, j);
1347
1348 return res;
1349 }
1350
1351
1352 // Matrix properties
1353
1354
1358 template<typename Matrix>
1359 inline bool is_square(const Matrix& m) {
1360 return (m.rows() == m.cols());
1361 }
1362
1363
1369 template<typename Matrix>
1371
1372 for (unsigned int i = 0; i < m.rows(); ++i)
1373 for (unsigned int j = 0; j < m.cols(); ++j)
1374 if(i != j && abs(m(i, j)) > tolerance)
1375 return false;
1376
1377 return true;
1378 }
1379
1385 template<typename Matrix>
1387
1388 if(!is_square(m))
1389 return false;
1390
1391 for (unsigned int i = 0; i < m.rows(); ++i)
1392 for (unsigned int j = 0; j < m.cols(); ++j)
1393 if(abs(m(i, j) - m(j, i)) > tolerance)
1394 return false;
1395
1396 return true;
1397 }
1398
1399
1405 template<typename Matrix>
1407
1408 if(!is_square(m))
1409 return false;
1410
1411 for (unsigned int i = 0; i < m.rows(); ++i)
1412 for (unsigned int j = i + 1; j < m.cols(); ++j)
1413 if (abs(m(i, j)) > tolerance)
1414 return false;
1415
1416 return true;
1417 }
1418
1419
1425 template<typename Matrix>
1427
1428 if(!is_square(m))
1429 return false;
1430
1431 for (unsigned int i = 0; i < m.rows(); ++i)
1432 for (unsigned int j = 0; j < i; ++j)
1433 if (abs(m(i, j)) > tolerance)
1434 return false;
1435
1436 return true;
1437 }
1438
1439
1449 template<typename Matrix, enable_matrix<Matrix> = true>
1451
1452 long unsigned int n = 0;
1453
1454 for (unsigned int i = 0; i < A.rows(); ++i)
1455 for (unsigned int j = 0; j < A.cols(); ++j)
1456 if (abs(A(i, j)) > tolerance)
1457 n++;
1458
1459 return (real(n) / A.rows()) / A.cols();
1460 }
1461
1462
1472 template<typename Matrix, enable_matrix<Matrix> = true>
1474
1475 return 1.0 - density(A, tolerance);
1476 }
1477
1478
1479 // Matrix decompositions
1480
1481
1489 template<typename Matrix1, typename Matrix2, typename Matrix3>
1490 inline void decompose_lu(const Matrix1& A, Matrix2& L, Matrix3& U) {
1491
1492 // Check the shapes of A, L and U
1493 unsigned int err = 0;
1494
1495 // Make sure to allocate L and U if they are empty
1496 L.resize(A.rows(), A.cols());
1497 U.resize(A.rows(), A.cols());
1498
1499 if (!is_square(A))
1500 err = A.rows();
1501
1502 if (A.rows() != L.rows())
1503 err = L.rows();
1504
1505 if (A.cols() != L.cols())
1506 err = L.cols();
1507
1508 if (A.rows() != U.rows())
1509 err = U.rows();
1510
1511 if (A.cols() != U.cols())
1512 err = U.cols();
1513
1514 if (err) {
1515 TH_MATH_ERROR("algebra::decompose_lu", err, MathError::InvalidArgument);
1517 return;
1518 }
1519
1520 using Type = matrix_element_t<Matrix1>;
1521
1522 // Set the diagonal of L to 1.0
1523 for (unsigned int i = 0; i < A.rows(); ++i)
1524 L(i, i) = (Type) 1.0;
1525
1526 // Compute L and U
1527 for(unsigned int i = 0; i < A.rows(); ++i) {
1528
1529 for(unsigned int j = 0; j < A.rows(); ++j) {
1530
1531 U(i, j) = A(i, j);
1532
1533 for(unsigned int k = 0; k < i; ++k)
1534 U(i, j) -= pair_inner_product(L(i, k), U(k, j));
1535 }
1536
1537 for(unsigned int j = i + 1; j < A.rows(); ++j) {
1538
1539 L(j, i) = A(j, i);
1540
1541 for(unsigned int k = 0; k < i; ++k)
1542 L(j, i) -= pair_inner_product(L(j, k), U(k, i));
1543
1544 L(j, i) /= U(i, i);
1545 }
1546 }
1547 }
1548
1549
1558 template<typename Matrix>
1560
1561 if (!is_square(A)) {
1562 TH_MATH_ERROR("algebra::decompose_lu_inplace", A.rows(), MathError::InvalidArgument);
1563 return mat_error(A);
1564 }
1565
1566 for(unsigned int j = 0; j < A.cols(); ++j) {
1567
1568 for(unsigned int i = j + 1; i < A.rows(); ++i) {
1569
1570 A(i, j) /= A(j, j);
1571
1572 for(unsigned int k = j + 1; k < A.rows(); ++k)
1573 A(i, k) -= pair_inner_product(A(i, j), A(j, k));
1574 }
1575 }
1576
1577 return A;
1578 }
1579
1580
1587 template<typename Matrix>
1589
1590 Matrix L;
1591 L.resize(A.rows(), A.cols());
1592 using Type = matrix_element_t<Matrix>;
1593
1594 if (!is_square(A)) {
1595 TH_MATH_ERROR("algebra::decompose_cholesky", A.rows(), MathError::InvalidArgument);
1596 return mat_error(L);
1597 }
1598
1599 if (!is_symmetric(A)) {
1600 TH_MATH_ERROR("algebra::decompose_cholesky", false, MathError::InvalidArgument);
1601 return mat_error(L);
1602 }
1603
1604 mat_zeroes(L);
1605
1606 for (unsigned int i = 0; i < A.cols(); ++i) {
1607
1608 for (unsigned int j = 0; j <= i; ++j) {
1609
1610 Type sum = pair_inner_product(L(i, 0), L(j, 0));
1611
1612 for (unsigned int k = 1; k < j; ++k)
1613 sum += pair_inner_product(L(i, k), L(j, k));
1614
1615 if (i == j) {
1616
1617 const Type sqr_diag = A(j, j) - sum;
1618
1619 // Additional check to ensure that the matrix is positive definite
1620 if (sqr_diag < MACH_EPSILON) {
1621 TH_MATH_ERROR("algebra::decompose_cholesky", sqr_diag, MathError::InvalidArgument);
1622 return mat_error(L);
1623 }
1624
1625 L(i, j) = sqrt(sqr_diag);
1626
1627 } else {
1628
1629 L(i, j) = (A(i, j) - sum) / L(j, j);
1630 }
1631 }
1632 }
1633
1634 return L;
1635 }
1636
1637
1643 template<typename Matrix>
1645
1646 using Type = matrix_element_t<Matrix>;
1647
1648 if (!is_square(A)) {
1649 TH_MATH_ERROR("algebra::decompose_cholesky_inplace", A.rows(), MathError::InvalidArgument);
1650 return mat_error(A);
1651 }
1652
1653 if (!is_symmetric(A)) {
1654 TH_MATH_ERROR("algebra::decompose_cholesky_inplace", false, MathError::InvalidArgument);
1655 return mat_error(A);
1656 }
1657
1658 // Compute the Cholesky decomposition in-place
1659 for (unsigned int k = 0; k < A.rows(); ++k) {
1660
1661 A(k, k) = sqrt(A(k, k));
1662
1663 for (unsigned int i = k + 1; i < A.rows(); ++i)
1664 A(i, k) = A(i, k) / A(k, k);
1665
1666 for (unsigned int j = k + 1; j < A.rows(); ++j)
1667 for (unsigned int i = j; i < A.cols(); ++i)
1668 A(i, j) -= pair_inner_product(A(i, k), A(j, k));
1669 }
1670
1671 // Zero out elements over the diagonal
1672 for (unsigned int i = 0; i < A.rows(); ++i)
1673 for (unsigned int j = i + 1; j < A.rows(); ++j)
1674 A(i, j) = Type(0.0);
1675
1676
1677 return A;
1678 }
1679
1680
1681 // Linear system solvers
1682
1683
1689 template<typename Matrix, typename Vector>
1690 inline Vector solve_triangular_lower(const Matrix& L, const Vector& b) {
1691
1692 Vector x;
1693 x.resize(L.cols());
1694 using Type = matrix_element_t<Matrix>;
1695
1696 if (!is_square(L)) {
1697 TH_MATH_ERROR("algebra::solve_triangular_lower", false, MathError::InvalidArgument);
1698 return vec_error(x);
1699 }
1700
1701 if (b.size() != L.rows()) {
1702 TH_MATH_ERROR("algebra::solve_triangular_lower", b.size(), MathError::InvalidArgument);
1703 return vec_error(x);
1704 }
1705
1706 // Solve using forward substitution
1707 for (unsigned int i = 0; i < L.cols(); ++i) {
1708
1709 Type sum = L(i, 0) * x[0];
1710
1711 for (unsigned int j = 1; j < i; ++j)
1712 sum += L(i, j) * x[j];
1713
1714 if (abs(L(i, i)) < MACH_EPSILON) {
1715 TH_MATH_ERROR("algebra::solve_triangular_lower", L(i, i), MathError::DivByZero);
1716 return vec_error(x);
1717 }
1718
1719 x[i] = (b[i] - sum) / L(i, i);
1720 }
1721
1722 return x;
1723 }
1724
1725
1731 template<typename Matrix, typename Vector>
1732 inline Vector solve_triangular_upper(const Matrix& U, const Vector& b) {
1733
1734 Vector x;
1735 x.resize(U.cols());
1736 using Type = matrix_element_t<Matrix>;
1737
1738 if (!is_square(U)) {
1739 TH_MATH_ERROR("solve_triangular_upper", false, MathError::InvalidArgument);
1740 return vec_error(x);
1741 }
1742
1743 if (b.size() != U.rows()) {
1744 TH_MATH_ERROR("solve_triangular_upper", b.size(), MathError::InvalidArgument);
1745 return vec_error(x);
1746 }
1747
1748 if (U.rows() == 1) {
1749 x[0] = b[0] / U(0, 0);
1750 return x;
1751 }
1752
1753 // Solve using backward substitution
1754 for (int i = U.rows() - 1; i >= 0; --i) {
1755
1756 Type sum = U(i, i + 1) * x[i + 1];
1757
1758 for (unsigned int j = i + 2; j < U.cols(); ++j)
1759 sum += U(i, j) * x[j];
1760
1761 if (abs(U(i, i)) < MACH_EPSILON) {
1762 TH_MATH_ERROR("solve_triangular_upper", U(i, i), MathError::DivByZero);
1763 return vec_error(x);
1764 }
1765
1766 x[i] = (b[i] - sum) / U(i, i);
1767 }
1768
1769 return x;
1770 }
1771
1772
1780 template<typename Matrix, typename Vector>
1781 inline Vector solve_triangular(const Matrix& T, const Vector& b) {
1782
1783 // Pick the correct solver
1785 return solve_triangular_lower(T, b);
1786 else if(is_upper_triangular(T))
1787 return solve_triangular_upper(T, b);
1788 else {
1789 Vector err;
1790 err.resize(b.size());
1791 return vec_error(err);
1792 }
1793 }
1794
1795
1806 template<typename Matrix, typename Vector>
1807 inline Vector& solve_lu_inplace(const Matrix& A, Vector& b) {
1808
1809 if (!is_square(A)) {
1810 TH_MATH_ERROR("algebra::solve_lu_inplace", A.rows(), MathError::InvalidArgument);
1811 return vec_error(b);
1812 }
1813
1814 if (A.rows() != b.size()) {
1815 TH_MATH_ERROR("algebra::solve_lu_inplace", A.rows(), MathError::InvalidArgument);
1816 return vec_error(b);
1817 }
1818
1819 using Type = matrix_element_t<Matrix>;
1820
1821 // Forward elimination: solves the lower system (L matrix).
1822 for (unsigned int i = 1; i < A.rows(); ++i) {
1823
1824 Type sum = Type(0.0);
1825
1826 for (unsigned int j = 0; j < i; ++j)
1827 sum += A(i, j) * b[j];
1828
1829 b[i] -= sum;
1830 }
1831
1832 // Backward elimination: solves the upper system (U matrix).
1833 for (int i = A.rows() - 1; i >= 0; --i) {
1834
1835 Type sum = Type(0.0);
1836
1837 for (unsigned int j = i + 1; j < A.rows(); ++j)
1838 sum += A(i, j) * b[j];
1839
1840 b[i] = (b[i] - sum) / A(i, i);
1841 }
1842
1843 return b;
1844 }
1845
1846
1854 template<typename Matrix, typename Vector>
1856
1857 // Apply in-place LU decomposition
1859
1860 // Apply forward and backward substitution
1861 return solve_lu_inplace(A, b);
1862 }
1863
1864
1875 template<typename Matrix1, typename Matrix2, typename Vector>
1876 inline Vector solve_lu(const Matrix1& L, const Matrix2& U, const Vector& b) {
1877
1878 Vector x = b;
1879
1880 if (!is_square(L)) {
1881 TH_MATH_ERROR("algebra::solve_lu", L.rows(), MathError::InvalidArgument);
1882 return vec_error(x);
1883 }
1884
1885 if (!is_square(U)) {
1886 TH_MATH_ERROR("algebra::solve_lu", U.rows(), MathError::InvalidArgument);
1887 return vec_error(x);
1888 }
1889
1890 if (L.rows() != U.rows()) {
1891 TH_MATH_ERROR("algebra::solve_lu", U.rows(), MathError::InvalidArgument);
1892 return vec_error(x);
1893 }
1894
1895 if (b.size() != L.rows()) {
1896 TH_MATH_ERROR("algebra::solve_lu", b.size(), MathError::InvalidArgument);
1897 return vec_error(x);
1898 }
1899
1900 using Type1 = matrix_element_t<Matrix1>;
1901
1902 // Forward elimination for L
1903 for (unsigned int i = 1; i < L.rows(); ++i) {
1904
1905 Type1 sum = Type1(0.0);
1906
1907 for (unsigned int j = 0; j < i; ++j)
1908 sum += L(i, j) * x[j];
1909
1910 x[i] -= sum;
1911 }
1912
1913 using Type2 = matrix_element_t<Matrix2>;
1914
1915 // Backward elimination for U
1916 for (int i = U.rows() - 1; i >= 0; --i) {
1917
1918 Type2 sum = Type2(0.0);
1919
1920 for (unsigned int j = i + 1; j < U.rows(); ++j)
1921 sum += U(i, j) * x[j];
1922
1923 if (abs(U(i, i)) < MACH_EPSILON) {
1924 TH_MATH_ERROR("algebra::solve_lu", U(i, i), MathError::DivByZero);
1925 return vec_error(x);
1926 }
1927
1928 x[i] = (x[i] - sum) / U(i, i);
1929 }
1930
1931 return x;
1932 }
1933
1934
1946 template<typename Matrix, typename Vector>
1947 inline Vector solve_cholesky(const Matrix& L, const Vector& b) {
1948
1949 Vector x = b;
1950
1951 if (!is_square(L)) {
1952 TH_MATH_ERROR("algebra::solve_cholesky", L.rows(), MathError::InvalidArgument);
1953 return vec_error(x);
1954 }
1955
1956 if (L.rows() != b.size()) {
1957 TH_MATH_ERROR("algebra::solve_cholesky", b.size(), MathError::InvalidArgument);
1958 return vec_error(x);
1959 }
1960
1961 using Type = matrix_element_t<Matrix>;
1962
1963 // Forward elimination for L
1964 for (unsigned int i = 0; i < L.rows(); ++i) {
1965
1966 Type sum = Type(0.0);
1967
1968 for (unsigned int j = 0; j < i; ++j)
1969 sum += L(i, j) * x[j];
1970
1971 x[i] = (x[i] - sum) / L(i, i);
1972 }
1973
1974 // Backward elimination for L transpose
1975 for (int i = L.rows() - 1; i >= 0; --i) {
1976
1977 Type sum = Type(0.0);
1978
1979 for (unsigned int j = i + 1; j < L.rows(); ++j)
1980 sum += L(j, i) * x[j];
1981
1982 x[i] = (x[i] - sum) / L(i, i);
1983 }
1984
1985 return x;
1986 }
1987
1988
1995 template<typename Matrix, typename Vector>
1996 inline Vector solve(const Matrix& A, const Vector& b) {
1997
1998 // Use LU decomposition to solve the linear system
1999 return solve_lu(A, b);
2000 }
2001
2002
2003 // Other composite operations
2004
2005
2011 template<typename Matrix1, typename Matrix2>
2013
2014 if(src.rows() != src.cols()) {
2015 TH_MATH_ERROR("algebra::inverse", src.rows(), MathError::InvalidArgument);
2016 return mat_error(dest);
2017 }
2018
2019 if(dest.rows() != src.rows()) {
2020 TH_MATH_ERROR("algebra::inverse", dest.rows(), MathError::InvalidArgument);
2021 return mat_error(dest);
2022 }
2023
2024 if(dest.cols() != src.cols()) {
2025 TH_MATH_ERROR("algebra::inverse", dest.cols(), MathError::InvalidArgument);
2026 return mat_error(dest);
2027 }
2028
2029 using Type = matrix_element_t<Matrix2>;
2030
2031 // Prepare extended matrix (A|B)
2032 Matrix1 A;
2033 A.resize(src.rows(), src.cols());
2034 dest.resize(src.rows(), src.cols());
2035 mat_copy(A, src);
2037
2038 // Iterate on all columns
2039 for (unsigned int i = 0; i < src.rows(); ++i) {
2040
2041 // Make sure the element on the diagonal
2042 // is non-zero by adding the first non-zero row
2043 if(A(i, i) == (Type) 0) {
2044
2045 bool flag = false;
2046
2047 // Iterate on all rows
2048 for (unsigned int j = i + 1; j < src.rows(); ++j) {
2049
2050 // Add the j-th row to the i-th row
2051 // if Aji is non-zero
2052 if(A(j, i) != (Type) 0) {
2053
2054 for (unsigned int k = 0; k < src.rows(); ++k) {
2055 A(i, k) += A(j, k);
2056 dest(i, k) += dest(j, k);
2057 }
2058
2059 flag = true;
2060 break;
2061 }
2062 }
2063
2064 if(!flag) {
2066 return mat_error(dest);
2067 }
2068 }
2069
2070 auto inv_pivot = ((Type) 1.0) / A(i, i);
2071
2072 // Use the current row to make all other
2073 // elements of the column equal to zero
2074 for (unsigned int j = 0; j < src.rows(); ++j) {
2075
2076 // Skip the current row
2077 if(j == i)
2078 continue;
2079
2080 // Multiplication coefficient for
2081 // the elision of Ajk
2082 const auto coeff = A(j, i) * inv_pivot;
2083
2084 for (unsigned int k = 0; k < src.rows(); ++k) {
2085 A(j, k) -= coeff * A(i, k);
2086 dest(j, k) -= coeff * dest(i, k);
2087 }
2088 }
2089
2090 // Divide the current row by the pivot
2091 for (unsigned int j = 0; j < src.rows(); ++j) {
2092 A(i, j) *= inv_pivot;
2093 dest(i, j) *= inv_pivot;
2094 }
2095
2096 }
2097
2098 return dest;
2099 }
2100
2101
2107 template<typename Matrix, typename MatrixInv = Matrix>
2108 inline MatrixInv inverse(const Matrix& m) {
2109 MatrixInv res;
2110 res.resize(m.rows(), m.cols());
2111 inverse(res, m);
2112 return res;
2113 }
2114
2115
2121 template<typename Matrix>
2122 inline Matrix& invert(Matrix& m) {
2123
2124 if(m.rows() != m.cols()) {
2125 TH_MATH_ERROR("algebra::invert", m.rows(), MathError::InvalidArgument);
2126 return mat_error(m);
2127 }
2128
2129 // Prepare extended matrix (A|B)
2130 Matrix tmp;
2131 tmp.resize(m.rows(), m.cols());
2132 inverse(tmp, m);
2133
2134 // Modify the matrix only when the inversion
2135 // has succeeded
2136 mat_copy(m, tmp);
2137 return m;
2138 }
2139
2140
2147 template<typename Matrix>
2148 inline auto det(const Matrix& A) {
2149
2150 Matrix LU = A;
2152
2153 // The determinant of a triangular matrix
2154 // is the product of the elements on its diagonal
2155 return diagonal_product(LU);
2156 }
2157
2158
2159 // Eigensolvers
2160
2161
2169 template<typename Matrix, typename Vector>
2170 inline auto rayleigh_quotient(const Matrix& A, const Vector& x) {
2171
2172 const auto p = dot(x, x);
2173
2174 // Check for division by zero
2175 if (abs(p) < MACH_EPSILON) {
2176 TH_MATH_ERROR("rayleigh_quotient", abs(p), MathError::DivByZero);
2177 return vector_element_t<Vector>(nan());
2178 }
2179
2180 return dot(x, transform(A, x)) / p;
2181 }
2182
2183
2195 template<typename Matrix, typename Vector>
2196 inline auto eigenvalue_power(
2197 const Matrix& A, const Vector& x,
2199
2200 using Type = matrix_element_t<Matrix>;
2201
2202 if (!is_square(A)) {
2203 TH_MATH_ERROR("eigenvalue_power", is_square(A), MathError::InvalidArgument);
2204 return make_error<Type>();
2205 }
2206
2207 if (x.size() != A.rows()) {
2208 TH_MATH_ERROR("eigenvalue_power", x.size(), MathError::InvalidArgument);
2209 return make_error<Type>();
2210 }
2211
2212 // Apply a first iteration to initialize
2213 // the current and previous vectors
2214 Vector x_prev = x;
2216
2217 // Iteration counter
2218 unsigned int i;
2219
2220 // Iteratively apply the matrix to the vector
2221 for (i = 1; i <= max_iter; ++i) {
2222
2223 x_prev = x_curr;
2225
2226 // Stop the algorithm when |x_k+1 +- x_k| is
2227 // less then the tolerance in module
2228 if (norm(x_curr - x_prev) <= tolerance ||
2230 break;
2231 }
2232
2233 // The algorithm did not converge
2234 if (i > max_iter) {
2235 TH_MATH_ERROR("eigenvalue_power", i, MathError::NoConvergence);
2236 return make_error<Type>();
2237 }
2238
2239 return dot(x_curr, A * x_curr);
2240 }
2241
2242
2256 template<typename Matrix, typename Vector1, typename Vector2 = Vector1>
2257 inline auto eigenpair_power(
2258 const Matrix& A, const Vector1& x, Vector2& v,
2260
2261 using Type = matrix_element_t<Matrix>;
2262
2263 if (!is_square(A)) {
2264 TH_MATH_ERROR("eigenpair_power", is_square(A), MathError::InvalidArgument);
2265 return make_error<Type>();
2266 }
2267
2268 if (x.size() != A.rows()) {
2269 TH_MATH_ERROR("eigenpair_power", x.size(), MathError::InvalidArgument);
2270 return make_error<Type>();
2271 }
2272
2273 if (v.size() != x.size()) {
2274 TH_MATH_ERROR("eigenpair_power", v.size(), MathError::InvalidArgument);
2275 return make_error<Type>();
2276 }
2277
2278 // Apply a first iteration to initialize
2279 // the current and previous vectors
2280 Vector1 x_prev = x;
2282
2283 // Iteration counter
2284 unsigned int i;
2285
2286 // Iteratively apply the matrix to the vector
2287 for (i = 1; i <= max_iter; ++i) {
2288
2289 x_prev = x_curr;
2291
2292 // Stop the algorithm when |x_k+1 +- x_k| is
2293 // less then the tolerance in module
2294 if (norm(x_curr - x_prev) <= tolerance ||
2296 break;
2297 }
2298
2299 // The algorithm did not converge
2300 if (i > max_iter) {
2301 TH_MATH_ERROR("eigenpair_power", i, MathError::NoConvergence);
2302 return make_error<Type>();
2303 }
2304
2305 // Overwrite with the eigenvector
2306 vec_copy(v, x_curr);
2307
2308 return dot(x_curr, A * x_curr);
2309 }
2310
2311
2323 template<typename Matrix, typename Vector>
2325 const Matrix& A, const Vector& x,
2327 unsigned int max_iter = ALGEBRA_EIGEN_ITER) {
2328
2329 using Type = matrix_element_t<Matrix>;
2330
2331 if (!is_square(A)) {
2332 TH_MATH_ERROR("eigenvalue_inverse", is_square(A), MathError::InvalidArgument);
2333 return make_error<Type>();
2334 }
2335
2336 if (A.rows() != x.size()) {
2337 TH_MATH_ERROR("eigenvalue_inverse", x.size(), MathError::InvalidArgument);
2338 return make_error<Type>();
2339 }
2340
2341 // Compute the LU decomposition of A to speed up system solution
2342 Matrix LU = A;
2344
2345 // Compute the first step to initialize the two vectors
2346 Vector x_prev = normalize(x);
2349
2350 // Iteration counter
2351 unsigned int i;
2352
2353 for (i = 1; i <= max_iter; ++i) {
2354
2355 // Solve the linear system using the LU decomposition
2356 x_prev = x_curr;
2358
2359 // Stop the algorithm when |x_k+1 +- x_k| is
2360 // less then the tolerance in module
2361 if (norm(x_curr - x_prev) <= tolerance ||
2363 break;
2364 }
2365
2366 // The algorithm did not converge
2367 if (i > max_iter) {
2368 TH_MATH_ERROR("eigenvalue_inverse", i, MathError::NoConvergence);
2369 return make_error<Type>();
2370 }
2371
2372 // A and inverse(A) have the same eigenvectors
2373 return dot(x_curr, A * x_curr);
2374 }
2375
2376
2390 template<typename Matrix, typename Vector1, typename Vector2>
2392 const Matrix& A, const Vector1& x, Vector2& v,
2394
2395 using Type = matrix_element_t<Matrix>;
2396
2397 if (!is_square(A)) {
2398 TH_MATH_ERROR("eigenpair_inverse", is_square(A), MathError::InvalidArgument);
2399 return make_error<Type>();
2400 }
2401
2402 if (A.rows() != x.size()) {
2403 TH_MATH_ERROR("eigenpair_inverse", x.size(), MathError::InvalidArgument);
2404 return make_error<Type>();
2405 }
2406
2407 // Compute the LU decomposition of A to speed up system solution
2408 Matrix LU = A;
2410
2411 // Compute the first step to initialize the two vectors
2415
2416 // Iteration counter
2417 unsigned int i;
2418
2419 for (i = 1; i <= max_iter; ++i) {
2420
2421 // Solve the linear system using the LU decomposition
2422 x_prev = x_curr;
2424
2425 // Stop the algorithm when |x_k+1 +- x_k| is
2426 // less then the tolerance in module
2427 if (norm(x_curr - x_prev) <= tolerance ||
2429 break;
2430 }
2431
2432 // The algorithm did not converge
2433 if (i > max_iter) {
2434 TH_MATH_ERROR("eigenpair_inverse", i, MathError::NoConvergence);
2435 return make_error<Type>();
2436 }
2437
2438 // Overwrite with the eigenvector
2439 vec_copy(v, x_curr);
2440
2441 // A and inverse(A) have the same eigenvectors
2442 return dot(x_curr, A * x_curr);
2443 }
2444
2445
2462 template<typename Matrix, typename Vector, typename T = matrix_element_t<Matrix>>
2464 const Matrix& A, const T& lambda, const Vector& x,
2466
2467 if (!is_square(A)) {
2468 TH_MATH_ERROR("eigenvector_inverse", is_square(A), MathError::InvalidArgument);
2469 return make_error<Vector>(A.rows());
2470 }
2471
2472 if (A.rows() != x.size()) {
2473 TH_MATH_ERROR("eigenvector_inverse", x.size(), MathError::InvalidArgument);
2474 return make_error<Vector>(A.rows());
2475 }
2476
2477 Matrix LU = A;
2479
2480 // Compute the LU decomposition of A
2481 // to speed up system solution
2483
2484 // Compute the first step to initialize the two vectors
2485 Vector v_prev = normalize(x);
2488
2489 // Iteration counter
2490 unsigned int i;
2491
2492 for (i = 1; i <= max_iter; ++i) {
2493
2494 // Solve the linear system using the LU decomposition
2495 v_prev = v_curr;
2497
2498 // Stop the algorithm when |x_k+1 +- x_k| is
2499 // less then the tolerance in module
2500 if (norm(v_curr - v_prev) <= tolerance ||
2502 break;
2503 }
2504
2505 // The algorithm did not converge
2506 if (i > max_iter) {
2507 TH_MATH_ERROR("eigenvector_inverse", i, MathError::NoConvergence);
2508 return vec_error(v_curr);
2509 }
2510
2511 return v_curr;
2512 }
2513
2514
2528 template<typename Matrix, typename Vector, typename T = matrix_element_t<Matrix>>
2530 const Matrix& A, const T& lambda, const Vector& x,
2532
2533 using Type = matrix_element_t<Matrix>;
2534
2535 if (!is_square(A)) {
2536 TH_MATH_ERROR("eigenvalue_rayleigh", is_square(A), MathError::InvalidArgument);
2537 return make_error<Type>();
2538 }
2539
2540 if (A.rows() != x.size()) {
2541 TH_MATH_ERROR("eigenvalue_rayleigh", x.size(), MathError::InvalidArgument);
2542 return make_error<Type>();
2543 }
2544
2545 // Keep track of the shifted matrix
2546 Matrix A_shift = A;
2548
2551
2552 Vector x_prev = normalize(x);
2554
2555 unsigned int i;
2556
2557 for (i = 1; i <= max_iter; ++i) {
2558
2559 // Update the eigenvalue approximation
2560 // using Rayleigh quotient (avoiding normalization)
2563
2564 // Shift the diagonal by the difference between
2565 // subsequent eigenvalues steps, to avoid copying matrix A
2567
2568 // Solve the linear system each time, as the shifted
2569 // matrix is continuously updated
2570 x_prev = x_curr;
2572
2573 // Stop the algorithm when |x_k+1 +- x_k| is
2574 // less then the tolerance in module
2575 if (norm(x_curr - x_prev) <= tolerance ||
2577 break;
2578 }
2579
2580 if (i > max_iter) {
2581 TH_MATH_ERROR("eigenvalue_rayleigh", i, MathError::NoConvergence);
2582 return make_error<Type>();
2583 }
2584
2585 return dot(x_curr, transform(A, x_curr));
2586 }
2587
2588
2604 template<typename Matrix, typename Vector1, typename Vector2,
2605 typename T = matrix_element_t<Matrix>>
2607 const Matrix& A, const T& lambda, const Vector1& x, Vector2& v,
2609
2610 using Type = matrix_element_t<Matrix>;
2611
2612 if (!is_square(A)) {
2613 TH_MATH_ERROR("eigenpair_rayleigh", is_square(A), MathError::InvalidArgument);
2614 return make_error<Type>();
2615 }
2616
2617 if (A.rows() != x.size()) {
2618 TH_MATH_ERROR("eigenpair_rayleigh", x.size(), MathError::InvalidArgument);
2619 return make_error<Type>();
2620 }
2621
2622 // Keep track of the shifted matrix
2623 Matrix A_shift = A;
2625
2628
2631
2632 unsigned int i;
2633
2634 for (i = 1; i <= max_iter; ++i) {
2635
2636 // Update the eigenvalue approximation
2637 // using Rayleigh quotient (avoiding normalization)
2640
2641 // Shift the diagonal by the difference between
2642 // subsequent eigenvalues steps, to avoid copying matrix A
2644
2645 // Solve the linear system each time, as the shifted
2646 // matrix is continuously updated
2647 x_prev = x_curr;
2649
2650 // Stop the algorithm when |x_k+1 +- x_k| is
2651 // less then the tolerance in module
2652 if (norm(x_curr - x_prev) <= tolerance ||
2654 break;
2655 }
2656
2657 if (i > max_iter) {
2658 TH_MATH_ERROR("eigenpair_rayleigh", i, MathError::NoConvergence);
2659 return make_error<Type>();
2660 }
2661
2662 // Overwrite with the eigenvector
2663 vec_copy(v, x_curr);
2664
2665 return dot(x_curr, transform(A, x_curr));
2666 }
2667 }
2668}
2669
2670#endif
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:219
bool is_lower_triangular(const Matrix &m, real tolerance=ALGEBRA_ELEMENT_TOL)
Returns whether the matrix is lower triangular.
Definition algebra.h:1406
Vector1 cross(const Vector1 &v1, const Vector2 &v2)
Compute the cross product between two tridimensional vectors.
Definition algebra.h:476
Matrix2 & mat_sum(Matrix1 &A, const Matrix2 &B)
Sum two matrices and store the result in the first matrix.
Definition algebra.h:858
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:1122
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:2257
Matrix & mat_scalmul(Field a, Matrix &m)
Multiply a matrix by a scalar of any compatible type.
Definition algebra.h:685
bool is_square(const Matrix &m)
Returns whether the matrix is square.
Definition algebra.h:1359
Matrix & invert(Matrix &m)
Invert the given matrix and overwrite it.
Definition algebra.h:2122
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:2196
auto trace(const Matrix &m)
Compute the trace of the given matrix.
Definition algebra.h:650
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:2391
Matrix & make_identity(Matrix &m)
Overwrite a matrix with the identity matrix.
Definition algebra.h:165
Matrix2 & mat_diff(Matrix1 &A, const Matrix2 &B)
Subtract two matrices and store the result in the first matrix.
Definition algebra.h:920
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:300
Vector & make_normalized(Vector &v)
Normalize a given vector overwriting it.
Definition algebra.h:430
Vector1 & vec_copy(Vector1 &dest, const Vector2 &src)
Copy a vector by overwriting another.
Definition algebra.h:241
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:1807
bool is_symmetric(const Matrix &m, real tolerance=ALGEBRA_ELEMENT_TOL)
Returns whether the matrix is symmetric.
Definition algebra.h:1386
Vector solve(const Matrix &A, const Vector &b)
Solve the linear system , finding using the best available algorithm.
Definition algebra.h:1996
ReturnVector vec_mat_mul(const Vector &v, const Matrix &A)
Multiply a vector by a matrix, returning the result.
Definition algebra.h:1298
Vector transform(const Matrix &A, const Vector &v)
Returns the matrix transformation of a vector.
Definition algebra.h:799
bool is_upper_triangular(const Matrix &m, real tolerance=ALGEBRA_ELEMENT_TOL)
Returns whether the matrix is upper triangular.
Definition algebra.h:1426
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:767
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:2463
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:1490
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:264
Vector solve_lu(Matrix A, Vector b)
Solve the linear system , finding .
Definition algebra.h:1855
Vector solve_triangular(const Matrix &T, const Vector &b)
Solve the linear system for triangular .
Definition algebra.h:1781
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:40
Matrix decompose_cholesky(const Matrix &A)
Decompose a symmetric positive definite matrix into a triangular matrix so that using Cholesky decom...
Definition algebra.h:1588
auto dot(const Vector1 &v, const Vector2 &w)
Computes the dot product between two vectors, using the Hermitian form if needed.
Definition algebra.h:454
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:1054
Matrix & mat_zeroes(Matrix &m)
Overwrite a matrix with all zeroes.
Definition algebra.h:181
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:355
auto diagonal_product(const Matrix &m)
Compute the product of the elements of the main diagonal of a generic matrix.
Definition algebra.h:668
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:2529
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:2606
auto rayleigh_quotient(const Matrix &A, const Vector &x)
Compute the Rayleigh quotient of a vector with respect to a matrix.
Definition algebra.h:2170
auto sqr_norm(const Vector &v)
Returns the square of the Euclidean/Hermitian norm of the given vector.
Definition algebra.h:378
auto det(const Matrix &A)
Compute the determinant of a square matrix.
Definition algebra.h:2148
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:336
Vector solve_triangular_lower(const Matrix &L, const Vector &b)
Solve the linear system for lower triangular .
Definition algebra.h:1690
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:983
Matrix & decompose_cholesky_inplace(Matrix &A)
Decompose a symmetric positive definite matrix in-place, overwriting the starting matrix,...
Definition algebra.h:1644
Vector2 & vec_sum(Vector1 &v1, const Vector2 &v2)
Sum two vectors and store the result in the first vector.
Definition algebra.h:1208
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:1154
Matrix1 & mat_copy(Matrix1 &dest, const Matrix2 &src)
Copy a matrix by overwriting another.
Definition algebra.h:213
Vector2 & vec_diff(Vector1 &v1, const Vector2 &v2)
Subtract two vectors and store the result in the first vector.
Definition algebra.h:1253
MatrixT hermitian(const Matrix &m)
Returns the hermitian of the given matrix.
Definition algebra.h:576
Matrix & make_transposed(Matrix &m)
Transpose the given matrix.
Definition algebra.h:521
Vector solve_triangular_upper(const Matrix &U, const Vector &b)
Solve the linear system for upper triangular .
Definition algebra.h:1732
Vector & vec_zeroes(Vector &v)
Overwrite a vector with all zeroes.
Definition algebra.h:197
Vector & vec_error(Vector &v)
Overwrite the given vector with the error vector with NaN values.
Definition algebra.h:58
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:1450
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:1947
Matrix1 & inverse(Matrix1 &dest, const Matrix2 &src)
Invert the given matrix.
Definition algebra.h:2012
auto norm(const Vector &v)
Returns the Euclidean/Hermitian norm of the given vector.
Definition algebra.h:395
bool is_diagonal(const Matrix &m, real tolerance=ALGEBRA_ELEMENT_TOL)
Returns whether the matrix is diagonal.
Definition algebra.h:1370
bool mat_equals(const Matrix1 &A, const Matrix2 &B, real tolerance=ALGEBRA_ELEMENT_TOL)
Checks whether two matrices are equal.
Definition algebra.h:1182
Vector normalize(const Vector &v)
Returns the normalized vector.
Definition algebra.h:405
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:2324
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:1473
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:1559
Matrix & make_hermitian(Matrix &m)
Compute the hermitian of a given matrix and overwrite it.
Definition algebra.h:595
Vector & vec_scalmul(Field a, Vector &v)
Multiply a vector by a scalar of any compatible type.
Definition algebra.h:727
MatrixT transpose(const Matrix &m)
Returns the transpose of the given matrix.
Definition algebra.h:504
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:347
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:90
dual2 abs(dual2 x)
Compute the absolute value of a second order dual number.
Definition dual2_functions.h:242
constexpr real ALGEBRA_ELEMENT_TOL
Tolerance for the elements of matrices.
Definition constants.h:279
Vector make_error()
Create a vector representing an error state, with all NaN values.
Definition algebra.h:103
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:215
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:74
@ 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