Theoretica
Mathematical Library
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
31 namespace algebra {
32
33
34 // Operations involving one matrix or vector
35
36
43 template<typename Matrix>
44 inline Matrix& mat_error(Matrix& m) {
45
46 using Type = matrix_element_t<Matrix>;
47
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);
51
52 return m;
53 }
54
55
61 template<typename Vector>
62 inline Vector& vec_error(Vector& v) {
63
64 using Type = vector_element_t<Vector>;
65
66 for (unsigned int i = 0; i < v.size(); ++i)
67 v[i] = Type(nan());
68
69 return v;
70 }
71
72
78 template<typename Vector>
79 inline Vector make_error(unsigned int n) {
80 Vector result;
81 result.resize(n);
82 return vec_error(result);
83 }
84
85
92 template<typename Matrix>
93 inline Matrix make_error(unsigned int rows, unsigned int cols) {
94 Matrix result;
95 result.resize(rows, cols);
96 return mat_error(result);
97 }
98
99
103 template<typename Matrix>
104 inline Matrix& make_identity(Matrix& m) {
105
106 using Type = matrix_element_t<Matrix>;
107
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);
111
112 return m;
113 }
114
115
119 template<typename Matrix>
120 inline Matrix& mat_zeroes(Matrix& m) {
121
122 using Type = matrix_element_t<Matrix>;
123
124 for (unsigned int i = 0; i < m.rows(); ++i)
125 for (unsigned int j = 0; j < m.cols(); ++j)
126 m(i, j) = (Type) 0;
127
128 return m;
129 }
130
131
135 template<typename Vector>
136 inline Vector& vec_zeroes(Vector& v) {
137
138 using Type = vector_element_t<Vector>;
139
140 for (unsigned int i = 0; i < v.size(); ++i)
141 v[i] = Type(0.0);
142
143 return v;
144 }
145
146
151 template<typename Matrix1, typename Matrix2>
152 inline Matrix1& mat_copy(Matrix1& dest, const Matrix2& src) {
153
154 dest.resize(src.rows(), src.cols());
155
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);
159
160 return dest;
161 }
162
163
169 template<typename Vector1, typename Vector2>
170 inline Vector1& vec_copy(Vector1& dest, const Vector2& src) {
171
172 dest.resize(src.size());
173
174 for (unsigned int i = 0; i < src.size(); ++i)
175 dest[i] = src[i];
176
177 return dest;
178 }
179
180
187 template<typename Matrix>
188 inline Matrix& mat_swap_rows(Matrix& A, unsigned int row1, unsigned int row2) {
189
190 using Type = matrix_element_t<Matrix>;
191
192 if (row1 >= A.rows()) {
193 TH_MATH_ERROR("algebra::mat_swap_rows", row1, MathError::InvalidArgument);
194 return mat_error(A);
195 }
196
197 if (row2 >= A.rows()) {
198 TH_MATH_ERROR("algebra::mat_swap_rows", row2, MathError::InvalidArgument);
199 return mat_error(A);
200 }
201
202 if (row1 == row2)
203 return A;
204
205 for (unsigned int j = 0; j < A.cols(); ++j) {
206
207 const Type tmp = A(row1, j);
208
209 A(row1, j) = A(row2, j);
210 A(row2, j) = tmp;
211 }
212
213 return A;
214 }
215
216
223 template<typename Matrix>
224 inline Matrix& mat_swap_cols(Matrix& A, unsigned int col1, unsigned int col2) {
225
226 using Type = matrix_element_t<Matrix>;
227
228 if (col1 >= A.cols()) {
229 TH_MATH_ERROR("algebra::mat_swap_cols", col1, MathError::InvalidArgument);
230 return mat_error(A);
231 }
232
233 if (col2 >= A.cols()) {
234 TH_MATH_ERROR("algebra::mat_swap_cols", col2, MathError::InvalidArgument);
235 return mat_error(A);
236 }
237
238 if (col1 == col2)
239 return A;
240
241 for (unsigned int i = 0; i < A.cols(); ++i) {
242
243 const Type tmp = A(i, col1);
244
245 A(i, col1) = A(i, col2);
246 A(i, col2) = tmp;
247 }
248
249 return A;
250 }
251
252
259 template<typename Matrix, typename Type = matrix_element_t<Matrix>>
260 inline Matrix& mat_shift_diagonal(Matrix& A, const Type& sigma) {
261
262 const unsigned int count = min(A.rows(), A.cols());
263
264 for (unsigned int i = 0; i < count; ++i)
265 A(i, i) += sigma;
266
267 return A;
268 }
269
270
278 template<typename Type>
279 inline auto pair_inner_product(const Type& v_i, const Type& w_i) {
280 return v_i * w_i;
281 }
282
290 template<typename Type>
291 inline auto pair_inner_product(const complex<Type>& v_i, const complex<Type>& w_i) {
292 return v_i * conjugate(w_i);
293 }
294
295
301 template<typename Vector>
302 inline auto sqr_norm(const Vector& v) {
303
304 auto sum = pair_inner_product(v[0], v[0]);
305
306 // Use conjugation for complex numbers
307 for (unsigned int i = 1; i < v.size(); ++i)
308 sum += pair_inner_product(v[i], v[i]);
309
310 return sum;
311 }
312
313
318 template<typename Vector>
319 inline auto norm(const Vector& v) {
320
321 return vector_element_t<Vector>(sqrt(sqr_norm(v)));
322 }
323
324
328 template<typename Vector>
329 inline Vector normalize(const Vector& v) {
330
331 Vector r;
332 r.resize(v.size());
333 vec_copy(r, v);
334
335 const auto m = norm(v);
336
337 if(abs(m) < MACH_EPSILON) {
338 TH_MATH_ERROR("algebra::normalize", m, MathError::DivByZero);
339 vec_error(r);
340 return r;
341 }
342
343 for (unsigned int i = 0; i < r.size(); ++i)
344 r[i] /= m;
345
346 return r;
347 }
348
349
353 template<typename Vector>
354 inline Vector& make_normalized(Vector& v) {
355
356 const auto m = norm(v);
357
358 if(abs(m) < MACH_EPSILON) {
359 TH_MATH_ERROR("algebra::make_normalized", m, MathError::DivByZero);
360 vec_error(v);
361 return v;
362 }
363
364 for (unsigned int i = 0; i < v.size(); ++i)
365 v[i] /= m;
366
367 return v;
368 }
369
370
377 template<typename Vector1, typename Vector2>
378 inline auto dot(const Vector1& v, const Vector2& w) {
379
380 if(v.size() != w.size()) {
381 TH_MATH_ERROR("algebra::dot", v.size(), MathError::InvalidArgument);
382 return vector_element_t<Vector1>(nan());
383 }
384
385 auto sum = pair_inner_product(v[0], w[0]);
386
387 // Use conjugation for complex numbers
388 for (unsigned int i = 1; i < v.size(); ++i)
389 sum += pair_inner_product(v[i], w[i]);
390
391 return sum;
392 }
393
394
399 template<typename Vector1, typename Vector2>
400 inline Vector1 cross(const Vector1& v1, const Vector2& v2) {
401
402 Vector1 v3;
403 v3.resize(3);
404
405 if(v1.size() != 3) {
406 TH_MATH_ERROR("algebra::cross", v1.size(), MathError::InvalidArgument);
407 return vec_error(v3);
408 }
409
410 if(v2.size() != 3) {
411 TH_MATH_ERROR("algebra::cross", v2.size(), MathError::InvalidArgument);
412 return vec_error(v3);
413 }
414
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];
418
419 return v3;
420 }
421
422
427 template<typename Matrix, typename MatrixT = Matrix>
428 inline MatrixT transpose(const Matrix& m) {
429
430 MatrixT res;
431 res.resize(m.cols(), m.rows());
432
433 for (unsigned int i = 0; i < m.rows(); ++i)
434 for (unsigned int j = 0; j < m.cols(); ++j)
435 res(j, i) = m(i, j);
436
437 return res;
438 }
439
440
444 template<typename Matrix>
445 inline Matrix& make_transposed(Matrix& m) {
446
447 if(m.rows() != m.cols()) {
448 TH_MATH_ERROR("algebra::make_transposed", m.rows(), MathError::InvalidArgument);
449 return mat_error(m);
450 }
451
452 for (unsigned int i = 0; i < m.rows(); ++i) {
453 for (unsigned int j = 0; j < i; ++j) {
454
455 const auto buff = m(i, j);
456 m(i, j) = m(j, i);
457 m(j, i) = buff;
458 }
459 }
460
461 return m;
462 }
463
464
471 template<typename Matrix1, typename Matrix2>
472 inline Matrix1& transpose(Matrix1& dest, const Matrix2& src) {
473
474 // Check that the two matrices have the correct
475 // number of rows and columns
476 if(src.rows() != dest.cols()) {
477 TH_MATH_ERROR("algebra::transpose", src.rows(), MathError::InvalidArgument);
478 return mat_error(dest);
479 }
480
481 if(src.cols() != dest.rows()) {
482 TH_MATH_ERROR("algebra::transpose", src.cols(), MathError::InvalidArgument);
483 return mat_error(dest);
484 }
485
486 // Overwrite dest with the transpose of src
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);
490
491 return dest;
492 }
493
494
499 template<typename Matrix, typename MatrixT = Matrix>
500 inline MatrixT hermitian(const Matrix& m) {
501
502 MatrixT res;
503 res.resize(m.cols(), m.rows());
504
505 for (unsigned int i = 0; i < m.rows(); ++i)
506 for (unsigned int j = 0; j < m.cols(); ++j)
507 res(j, i) = conjugate(m(i, j));
508
509 return res;
510 }
511
512
518 template<typename Matrix>
519 inline Matrix& make_hermitian(Matrix& m) {
520
521 if(m.rows() != m.cols()) {
522 TH_MATH_ERROR("algebra::hermitian", m.rows(), MathError::InvalidArgument);
523 return mat_error(m);
524 }
525
526 for (unsigned int i = 0; i < m.rows(); ++i) {
527 for (unsigned int j = 0; j < i; ++j) {
528
529 const auto buff = m(i, j);
530 m(i, j) = conjugate(m(j, i));
531 m(j, i) = conjugate(buff);
532 }
533 }
534
535 return m;
536 }
537
538
546 template<typename Matrix1, typename Matrix2>
547 inline Matrix1& hermitian(Matrix1& dest, const Matrix2& src) {
548
549 // Check that the two matrices have the correct
550 // number of rows and columns
551 if(src.rows() != dest.cols()) {
552 TH_MATH_ERROR("algebra::hermitian", src.rows(), MathError::InvalidArgument);
553 return mat_error(dest);
554 }
555
556 if(src.cols() != dest.rows()) {
557 TH_MATH_ERROR("algebra::hermitian", src.cols(), MathError::InvalidArgument);
558 return mat_error(dest);
559 }
560
561 // Overwrite dest with the transpose of src
562 for (unsigned int i = 0; i < src.rows(); ++i)
563 for (unsigned int j = 0; j < src.cols(); ++j)
564 dest(j, i) = conjugate(src(i, j));
565
566 return dest;
567 }
568
569
573 template<typename Matrix>
574 inline auto trace(const Matrix& m) {
575
576 auto sum = m(0, 0);
577 const size_t n = min(m.rows(), m.cols());
578
579 for (unsigned int i = 1; i < n; ++i)
580 sum += m(i, i);
581
582 return sum;
583 }
584
585
591 template<typename Matrix>
592 inline auto diagonal_product(const Matrix& m) {
593
594 auto mul = m(0, 0);
595 const size_t n = min(m.rows(), m.cols());
596
597 for (unsigned int i = 1; i < n; ++i)
598 mul *= m(i, i);
599
600 return mul;
601 }
602
603
608 template<typename Field, typename Matrix>
609 inline Matrix& mat_scalmul(Field a, Matrix& m) {
610
611 for (unsigned int i = 0; i < m.rows(); ++i)
612 for (unsigned int j = 0; j < m.cols(); ++j)
613 m(i, j) *= a;
614
615 return m;
616 }
617
618
625 template<typename Field, typename Matrix1, typename Matrix2>
626 inline Matrix1& mat_scalmul(Matrix1& dest, Field a, const Matrix2& src) {
627
628 if(src.rows() != dest.rows()) {
629 TH_MATH_ERROR("algebra::mat_scalmul", src.rows(), MathError::InvalidArgument);
630 return mat_error(dest);
631 }
632
633 if(src.cols() != dest.cols()) {
634 TH_MATH_ERROR("algebra::mat_scalmul", src.cols(), MathError::InvalidArgument);
635 return mat_error(dest);
636 }
637
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);
641
642 return dest;
643 }
644
645
650 template<typename Field, typename Vector>
651 inline Vector& vec_scalmul(Field a, Vector& v) {
652
653 for (unsigned int i = 0; i < v.size(); ++i)
654 v[i] *= a;
655
656 return v;
657 }
658
659
666 template<typename Field, typename Vector1, typename Vector2>
667 inline Vector1& vec_scalmul(Vector1& dest, Field a, const Vector2& src) {
668
669 if(src.size() != dest.size()) {
670 TH_MATH_ERROR("algebra::vec_scalmul", src.size(), MathError::InvalidArgument);
671 return vec_error(dest);
672 }
673
674 for (unsigned int i = 0; i < src.size(); ++i)
675 dest[i] = a * src[i];
676
677 return dest;
678 }
679
680
681 // Operations involving a matrix and a vector
682
683
690 template<typename Matrix, typename Vector>
691 inline Vector& apply_transform(const Matrix& A, Vector& v) {
692
693 if(v.size() != A.cols()) {
694 TH_MATH_ERROR("algebra::apply_transform", v.size(), MathError::InvalidArgument);
695 return vec_error(v);
696 }
697
698 Vector res;
699 res.resize(v.size());
700 vec_zeroes(res);
701
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];
705
706 vec_copy(res, v);
707 return v;
708 }
709
710
716 template<typename Matrix, typename Vector>
717 inline Vector transform(const Matrix& A, const Vector& v) {
718
719 Vector res;
720 res.resize(v.size());
721
722 if(v.size() != A.cols()) {
723 TH_MATH_ERROR("algebra::transform", v.size(), MathError::InvalidArgument);
724 return vec_error(res);
725 }
726
727 vec_zeroes(res);
728
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];
732
733 return res;
734 }
735
736
744 template<typename Matrix, typename Vector1, typename Vector2>
745 inline Vector1& transform(Vector1& res, const Matrix& A, const Vector2& v) {
746
747 if(v.size() != A.cols()) {
748 TH_MATH_ERROR("algebra::transform", v.size(), MathError::InvalidArgument);
749 return vec_error(res);
750 }
751
752 if(res.size() != v.size()) {
753 TH_MATH_ERROR("algebra::transform", res.size(), MathError::InvalidArgument);
754 return vec_error(res);
755 }
756
757 vec_zeroes(res);
758
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];
762
763 return res;
764 }
765
766
767 // Operations involving multiple matrices or vectors
768
769
775 template<typename Matrix1, typename Matrix2>
776 inline Matrix2& mat_sum(Matrix1& A, const Matrix2& B) {
777
778 if(A.rows() != B.rows()) {
779 TH_MATH_ERROR("algebra::mat_sum", A.rows(), MathError::InvalidArgument);
780 return mat_error(A);
781 }
782
783 if(A.cols() != B.cols()) {
784 TH_MATH_ERROR("algebra::mat_sum", A.cols(), MathError::InvalidArgument);
785 return mat_error(A);
786 }
787
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);
791
792 return A;
793 }
794
795
801 template<typename Matrix1, typename Matrix2, typename Matrix3>
802 inline Matrix1& mat_sum(Matrix1& res, const Matrix2& A, const Matrix3& B) {
803
804 if(A.rows() != B.rows()) {
805 TH_MATH_ERROR("algebra::mat_sum", A.rows(), MathError::InvalidArgument);
806 return mat_error(res);
807 }
808
809 if(A.cols() != B.cols()) {
810 TH_MATH_ERROR("algebra::mat_sum", A.cols(), MathError::InvalidArgument);
811 return mat_error(res);
812 }
813
814 if(res.rows() != A.rows()) {
815 TH_MATH_ERROR("algebra::mat_sum", res.rows(), MathError::InvalidArgument);
816 return mat_error(res);
817 }
818
819 if(res.cols() != A.cols()) {
820 TH_MATH_ERROR("algebra::mat_sum", res.cols(), MathError::InvalidArgument);
821 return mat_error(res);
822 }
823
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);
827
828 return res;
829 }
830
831
837 template<typename Matrix1, typename Matrix2>
838 inline Matrix2& mat_diff(Matrix1& A, const Matrix2& B) {
839
840 if(A.rows() != B.rows()) {
841 TH_MATH_ERROR("algebra::mat_diff", A.rows(), MathError::InvalidArgument);
842 return mat_error(A);
843 }
844
845 if(A.cols() != B.cols()) {
846 TH_MATH_ERROR("algebra::mat_diff", A.cols(), MathError::InvalidArgument);
847 return mat_error(A);
848 }
849
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);
853
854 return A;
855 }
856
857
863 template<typename Matrix1, typename Matrix2, typename Matrix3>
864 inline Matrix1& mat_diff(Matrix1& res, const Matrix2& A, const Matrix3& B) {
865
866 if(A.rows() != B.rows()) {
867 TH_MATH_ERROR("algebra::mat_diff", A.rows(), MathError::InvalidArgument);
868 return mat_error(res);
869 }
870
871 if(A.cols() != B.cols()) {
872 TH_MATH_ERROR("algebra::mat_diff", A.cols(), MathError::InvalidArgument);
873 return mat_error(res);
874 }
875
876 if(res.rows() != A.rows()) {
877 TH_MATH_ERROR("algebra::mat_diff", res.rows(), MathError::InvalidArgument);
878 return mat_error(res);
879 }
880
881 if(res.cols() != A.cols()) {
882 TH_MATH_ERROR("algebra::mat_diff", res.cols(), MathError::InvalidArgument);
883 return mat_error(res);
884 }
885
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);
889
890 return res;
891 }
892
893
900 template<typename Field1, typename Matrix1, typename Field2, typename Matrix2>
901 inline Matrix2& mat_lincomb(
902 Field1 alpha, Matrix1& A, Field2 beta, const Matrix2& B) {
903
904 if(A.rows() != B.rows()) {
905 TH_MATH_ERROR("algebra::mat_lincomb", A.rows(), MathError::InvalidArgument);
906 return mat_error(A);
907 }
908
909 if(A.cols() != B.cols()) {
910 TH_MATH_ERROR("algebra::mat_lincomb", A.cols(), MathError::InvalidArgument);
911 return mat_error(A);
912 }
913
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;
917
918 return A;
919 }
920
921
931 template<typename Matrix1, typename Field1, typename Matrix2,
932 typename Field2, typename Matrix3>
933 inline Matrix1& mat_lincomb(
934 Matrix1& res, Field1 alpha, const Matrix2& A, Field2 beta, const Matrix3& B) {
935
936 if(A.rows() != B.rows()) {
937 TH_MATH_ERROR("algebra::mat_lincomb", A.rows(), MathError::InvalidArgument);
938 return mat_error(res);
939 }
940
941 if(A.cols() != B.cols()) {
942 TH_MATH_ERROR("algebra::mat_lincomb", A.cols(), MathError::InvalidArgument);
943 return mat_error(res);
944 }
945
946 if(res.rows() != A.rows()) {
947 TH_MATH_ERROR("algebra::mat_lincomb", res.rows(), MathError::InvalidArgument);
948 return mat_error(res);
949 return res;
950 }
951
952 if(res.cols() != A.cols()) {
953 TH_MATH_ERROR("algebra::mat_lincomb", res.cols(), MathError::InvalidArgument);
954 return mat_error(res);
955 }
956
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;
960
961 return A;
962 }
963
964
971 template<typename Matrix1, typename Matrix2, typename Matrix3 = Matrix1>
972 inline Matrix3 mat_mul(const Matrix1& A, const Matrix2& B) {
973
974 Matrix3 R;
975 R.resize(A.rows(), B.cols());
976
977 if(A.cols() != B.rows()) {
978 TH_MATH_ERROR("algebra::mat_mul", A.cols(), MathError::InvalidArgument);
979 return mat_error(R);
980 }
981
982 mat_zeroes(R);
983
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);
988
989 return R;
990 }
991
992
1000 template<typename Matrix1, typename Matrix2, typename Matrix3>
1001 inline Matrix1& mat_mul(Matrix1& R, const Matrix2& A, const Matrix3& B) {
1002
1003 if(R.rows() != A.rows()) {
1004 TH_MATH_ERROR("algebra::mat_mul", R.rows(), MathError::InvalidArgument);
1005 return mat_error(R);
1006 }
1007
1008 if(R.cols() != B.cols()) {
1009 TH_MATH_ERROR("algebra::mat_mul", R.cols(), MathError::InvalidArgument);
1010 return mat_error(R);
1011 }
1012
1013 if(A.cols() != B.rows()) {
1014 TH_MATH_ERROR("algebra::mat_mul", A.cols(), MathError::InvalidArgument);
1015 return mat_error(R);
1016 }
1017
1018 mat_zeroes(R);
1019
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);
1024
1025 return R;
1026 }
1027
1028
1039 template<typename Matrix1, typename Matrix2, typename Matrix3 = Matrix1>
1040 inline Matrix3 mat_transpose_mul(const Matrix1& A, const Matrix2& B) {
1041
1042 Matrix3 R;
1043 R.resize(A.cols(), B.cols());
1044
1045 if(A.rows() != B.rows()) {
1046 TH_MATH_ERROR("algebra::mat_transpose_mul", A.rows(), MathError::InvalidArgument);
1047 return mat_error(R);
1048 }
1049
1050 mat_zeroes(R);
1051
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);
1056
1057 return R;
1058 }
1059
1060
1071 template<typename Matrix1, typename Matrix2, typename Matrix3 = Matrix1>
1072 inline Matrix3 mat_mul_transpose(const Matrix1& A, const Matrix2& B) {
1073
1074 Matrix3 R;
1075 R.resize(A.rows(), B.rows());
1076
1077 if(A.cols() != B.cols()) {
1078 TH_MATH_ERROR("algebra::mat_mul_transpose", A.rows(), MathError::InvalidArgument);
1079 return mat_error(R);
1080 }
1081
1082 mat_zeroes(R);
1083
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);
1088
1089 return R;
1090 }
1091
1092
1099 template<typename Matrix1, typename Matrix2>
1100 inline bool mat_equals(
1101 const Matrix1& A, const Matrix2& B, real tolerance = ALGEBRA_ELEMENT_TOL) {
1102
1103 if(A.rows() != B.rows() || A.cols() != B.cols())
1104 return false;
1105
1106 for (unsigned int i = 0; i < A.rows(); ++i) {
1107 for (unsigned int j = 0; j < A.cols(); ++j) {
1108
1109 const auto diff = abs(A(i, j) - B(i, j));
1110
1111 if(diff > tolerance || is_nan(diff))
1112 return false;
1113 }
1114 }
1115
1116 return true;
1117 }
1118
1119
1125 template<typename Vector1, typename Vector2>
1126 inline Vector2& vec_sum(Vector1& v1, const Vector2& v2) {
1127
1128 if(v1.size() != v2.size()) {
1129 TH_MATH_ERROR("algebra::vec_sum", v1.size(), MathError::InvalidArgument);
1130 return vec_error(v1);
1131 }
1132
1133 for (unsigned int i = 0; i < v1.size(); ++i)
1134 v1[i] = v1[i] + v2[i];
1135
1136 return v1;
1137 }
1138
1139
1145 template<typename Vector1, typename Vector2, typename Vector3>
1146 inline Vector1& vec_sum(Vector1& res, const Vector2& v1, const Vector3& v2) {
1147
1148 if(v1.size() != v2.size()) {
1149 TH_MATH_ERROR("algebra::vec_sum", v1.size(), MathError::InvalidArgument);
1150 return vec_error(res);
1151 }
1152
1153 if(res.size() != v1.size()) {
1154 TH_MATH_ERROR("algebra::vec_sum", res.size(), MathError::InvalidArgument);
1155 return vec_error(res);
1156 }
1157
1158 for (unsigned int i = 0; i < v1.size(); ++i)
1159 res[i] = v1[i] + v2[i];
1160
1161 return res;
1162 }
1163
1164
1170 template<typename Vector1, typename Vector2>
1171 inline Vector2& vec_diff(Vector1& v1, const Vector2& v2) {
1172
1173 if(v1.size() != v2.size()) {
1174 TH_MATH_ERROR("algebra::vec_diff", v1.size(), MathError::InvalidArgument);
1175 return vec_error(v1);
1176 }
1177
1178 for (unsigned int i = 0; i < v1.size(); ++i)
1179 v1[i] = v1[i] - v2[i];
1180
1181 return v1;
1182 }
1183
1184
1190 template<typename Vector1, typename Vector2, typename Vector3>
1191 inline Vector1& vec_diff(Vector1& res, const Vector2& v1, const Vector3& v2) {
1192
1193 if(v1.size() != v2.size()) {
1194 TH_MATH_ERROR("algebra::vec_diff", v1.size(), MathError::InvalidArgument);
1195 return vec_error(res);
1196 }
1197
1198 if(res.size() != v1.size()) {
1199 TH_MATH_ERROR("algebra::vec_diff", res.size(), MathError::InvalidArgument);
1200 return vec_error(res);
1201 }
1202
1203 for (unsigned int i = 0; i < v1.size(); ++i)
1204 res[i] = v1[i] - v2[i];
1205
1206 return res;
1207 }
1208
1209
1210 // Matrix properties
1211
1212
1216 template<typename Matrix>
1217 inline bool is_square(const Matrix& m) {
1218 return (m.rows() == m.cols());
1219 }
1220
1221
1227 template<typename Matrix>
1228 inline bool is_diagonal(const Matrix& m, real tolerance = ALGEBRA_ELEMENT_TOL) {
1229
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)
1233 return false;
1234
1235 return true;
1236 }
1237
1243 template<typename Matrix>
1244 inline bool is_symmetric(const Matrix& m, real tolerance = ALGEBRA_ELEMENT_TOL) {
1245
1246 if(!is_square(m))
1247 return false;
1248
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)
1252 return false;
1253
1254 return true;
1255 }
1256
1257
1263 template<typename Matrix>
1264 inline bool is_lower_triangular(const Matrix& m, real tolerance = ALGEBRA_ELEMENT_TOL) {
1265
1266 if(!is_square(m))
1267 return false;
1268
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)
1272 return false;
1273
1274 return true;
1275 }
1276
1277
1283 template<typename Matrix>
1284 inline bool is_upper_triangular(const Matrix& m, real tolerance = ALGEBRA_ELEMENT_TOL) {
1285
1286 if(!is_square(m))
1287 return false;
1288
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)
1292 return false;
1293
1294 return true;
1295 }
1296
1297
1307 template<typename Matrix, enable_matrix<Matrix> = true>
1308 inline real density(const Matrix& A, real tolerance = ALGEBRA_ELEMENT_TOL) {
1309
1310 long unsigned int n = 0;
1311
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)
1315 n++;
1316
1317 return (real(n) / A.rows()) / A.cols();
1318 }
1319
1320
1330 template<typename Matrix, enable_matrix<Matrix> = true>
1331 inline real sparsity(const Matrix& A, real tolerance = ALGEBRA_ELEMENT_TOL) {
1332
1333 return 1.0 - density(A, tolerance);
1334 }
1335
1336
1337 // Matrix decompositions
1338
1339
1347 template<typename Matrix1, typename Matrix2, typename Matrix3>
1348 inline void decompose_lu(const Matrix1& A, Matrix2& L, Matrix3& U) {
1349
1350 // Check the shapes of A, L and U
1351 unsigned int err = 0;
1352
1353 // Make sure to allocate L and U if they are empty
1354 L.resize(A.rows(), A.cols());
1355 U.resize(A.rows(), A.cols());
1356
1357 if (!is_square(A))
1358 err = A.rows();
1359
1360 if (A.rows() != L.rows())
1361 err = L.rows();
1362
1363 if (A.cols() != L.cols())
1364 err = L.cols();
1365
1366 if (A.rows() != U.rows())
1367 err = U.rows();
1368
1369 if (A.cols() != U.cols())
1370 err = U.cols();
1371
1372 if (err) {
1373 TH_MATH_ERROR("algebra::decompose_lu", err, MathError::InvalidArgument);
1374 mat_error(L); mat_error(U);
1375 return;
1376 }
1377
1378 using Type = matrix_element_t<Matrix1>;
1379
1380 // Set the diagonal of L to 1.0
1381 for (unsigned int i = 0; i < A.rows(); ++i)
1382 L(i, i) = (Type) 1.0;
1383
1384 // Compute L and U
1385 for(unsigned int i = 0; i < A.rows(); ++i) {
1386
1387 for(unsigned int j = 0; j < A.rows(); ++j) {
1388
1389 U(i, j) = A(i, j);
1390
1391 for(unsigned int k = 0; k < i; ++k)
1392 U(i, j) -= pair_inner_product(L(i, k), U(k, j));
1393 }
1394
1395 for(unsigned int j = i + 1; j < A.rows(); ++j) {
1396
1397 L(j, i) = A(j, i);
1398
1399 for(unsigned int k = 0; k < i; ++k)
1400 L(j, i) -= pair_inner_product(L(j, k), U(k, i));
1401
1402 L(j, i) /= U(i, i);
1403 }
1404 }
1405 }
1406
1407
1416 template<typename Matrix>
1417 inline Matrix& decompose_lu_inplace(Matrix& A) {
1418
1419 if (!is_square(A)) {
1420 TH_MATH_ERROR("algebra::decompose_lu_inplace", A.rows(), MathError::InvalidArgument);
1421 return mat_error(A);
1422 }
1423
1424 for(unsigned int j = 0; j < A.cols(); ++j) {
1425
1426 for(unsigned int i = j + 1; i < A.rows(); ++i) {
1427
1428 A(i, j) /= A(j, j);
1429
1430 for(unsigned int k = j + 1; k < A.rows(); ++k)
1431 A(i, k) -= pair_inner_product(A(i, j), A(j, k));
1432 }
1433 }
1434
1435 return A;
1436 }
1437
1438
1445 template<typename Matrix>
1446 inline Matrix decompose_cholesky(const Matrix& A) {
1447
1448 Matrix L;
1449 L.resize(A.rows(), A.cols());
1450 using Type = matrix_element_t<Matrix>;
1451
1452 if (!is_square(A)) {
1453 TH_MATH_ERROR("algebra::decompose_cholesky", A.rows(), MathError::InvalidArgument);
1454 return mat_error(L);
1455 }
1456
1457 if (!is_symmetric(A)) {
1458 TH_MATH_ERROR("algebra::decompose_cholesky", false, MathError::InvalidArgument);
1459 return mat_error(L);
1460 }
1461
1462 mat_zeroes(L);
1463
1464 for (unsigned int i = 0; i < A.cols(); ++i) {
1465
1466 for (unsigned int j = 0; j <= i; ++j) {
1467
1468 Type sum = pair_inner_product(L(i, 0), L(j, 0));
1469
1470 for (unsigned int k = 1; k < j; ++k)
1471 sum += pair_inner_product(L(i, k), L(j, k));
1472
1473 if (i == j) {
1474
1475 const Type sqr_diag = A(j, j) - sum;
1476
1477 // Additional check to ensure that the matrix is positive definite
1478 if (sqr_diag < MACH_EPSILON) {
1479 TH_MATH_ERROR("algebra::decompose_cholesky", sqr_diag, MathError::InvalidArgument);
1480 return mat_error(L);
1481 }
1482
1483 L(i, j) = sqrt(sqr_diag);
1484
1485 } else {
1486
1487 L(i, j) = (A(i, j) - sum) / L(j, j);
1488 }
1489 }
1490 }
1491
1492 return L;
1493 }
1494
1495
1501 template<typename Matrix>
1502 inline Matrix& decompose_cholesky_inplace(Matrix& A) {
1503
1504 using Type = matrix_element_t<Matrix>;
1505
1506 if (!is_square(A)) {
1507 TH_MATH_ERROR("algebra::decompose_cholesky_inplace", A.rows(), MathError::InvalidArgument);
1508 return mat_error(A);
1509 }
1510
1511 if (!is_symmetric(A)) {
1512 TH_MATH_ERROR("algebra::decompose_cholesky_inplace", false, MathError::InvalidArgument);
1513 return mat_error(A);
1514 }
1515
1516 // Compute the Cholesky decomposition in-place
1517 for (unsigned int k = 0; k < A.rows(); ++k) {
1518
1519 A(k, k) = sqrt(A(k, k));
1520
1521 for (unsigned int i = k + 1; i < A.rows(); ++i)
1522 A(i, k) = A(i, k) / A(k, k);
1523
1524 for (unsigned int j = k + 1; j < A.rows(); ++j)
1525 for (unsigned int i = j; i < A.cols(); ++i)
1526 A(i, j) -= pair_inner_product(A(i, k), A(j, k));
1527 }
1528
1529 // Zero out elements over the diagonal
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);
1533
1534
1535 return A;
1536 }
1537
1538
1539 // Linear system solvers
1540
1541
1547 template<typename Matrix, typename Vector>
1548 inline Vector solve_triangular_lower(const Matrix& L, const Vector& b) {
1549
1550 Vector x;
1551 x.resize(L.cols());
1552 using Type = matrix_element_t<Matrix>;
1553
1554 if (!is_square(L)) {
1555 TH_MATH_ERROR("algebra::solve_triangular_lower", false, MathError::InvalidArgument);
1556 return vec_error(x);
1557 }
1558
1559 if (b.size() != L.rows()) {
1560 TH_MATH_ERROR("algebra::solve_triangular_lower", b.size(), MathError::InvalidArgument);
1561 return vec_error(x);
1562 }
1563
1564 // Solve using forward substitution
1565 for (unsigned int i = 0; i < L.cols(); ++i) {
1566
1567 Type sum = L(i, 0) * x[0];
1568
1569 for (unsigned int j = 1; j < i; ++j)
1570 sum += L(i, j) * x[j];
1571
1572 if (abs(L(i, i)) < MACH_EPSILON) {
1573 TH_MATH_ERROR("algebra::solve_triangular_lower", L(i, i), MathError::DivByZero);
1574 return vec_error(x);
1575 }
1576
1577 x[i] = (b[i] - sum) / L(i, i);
1578 }
1579
1580 return x;
1581 }
1582
1583
1589 template<typename Matrix, typename Vector>
1590 inline Vector solve_triangular_upper(const Matrix& U, const Vector& b) {
1591
1592 Vector x;
1593 x.resize(U.cols());
1594 using Type = matrix_element_t<Matrix>;
1595
1596 if (!is_square(U)) {
1597 TH_MATH_ERROR("solve_triangular_upper", false, MathError::InvalidArgument);
1598 return vec_error(x);
1599 }
1600
1601 if (b.size() != U.rows()) {
1602 TH_MATH_ERROR("solve_triangular_upper", b.size(), MathError::InvalidArgument);
1603 return vec_error(x);
1604 }
1605
1606 if (U.rows() == 1) {
1607 x[0] = b[0] / U(0, 0);
1608 return x;
1609 }
1610
1611 // Solve using backward substitution
1612 for (int i = U.rows() - 1; i >= 0; --i) {
1613
1614 Type sum = U(i, i + 1) * x[i + 1];
1615
1616 for (unsigned int j = i + 2; j < U.cols(); ++j)
1617 sum += U(i, j) * x[j];
1618
1619 if (abs(U(i, i)) < MACH_EPSILON) {
1620 TH_MATH_ERROR("solve_triangular_upper", U(i, i), MathError::DivByZero);
1621 return vec_error(x);
1622 }
1623
1624 x[i] = (b[i] - sum) / U(i, i);
1625 }
1626
1627 return x;
1628 }
1629
1630
1638 template<typename Matrix, typename Vector>
1639 inline Vector solve_triangular(const Matrix& T, const Vector& b) {
1640
1641 // Pick the correct solver
1642 if (is_lower_triangular(T))
1643 return solve_triangular_lower(T, b);
1644 else if(is_upper_triangular(T))
1645 return solve_triangular_upper(T, b);
1646 else {
1647 Vector err;
1648 err.resize(b.size());
1649 return vec_error(err);
1650 }
1651 }
1652
1653
1664 template<typename Matrix, typename Vector>
1665 inline Vector& solve_lu_inplace(const Matrix& A, Vector& b) {
1666
1667 if (!is_square(A)) {
1668 TH_MATH_ERROR("algebra::solve_lu_inplace", A.rows(), MathError::InvalidArgument);
1669 return vec_error(b);
1670 }
1671
1672 if (A.rows() != b.size()) {
1673 TH_MATH_ERROR("algebra::solve_lu_inplace", A.rows(), MathError::InvalidArgument);
1674 return vec_error(b);
1675 }
1676
1677 using Type = matrix_element_t<Matrix>;
1678
1679 // Forward elimination: solves the lower system (L matrix).
1680 for (unsigned int i = 1; i < A.rows(); ++i) {
1681
1682 Type sum = Type(0.0);
1683
1684 for (unsigned int j = 0; j < i; ++j)
1685 sum += A(i, j) * b[j];
1686
1687 b[i] -= sum;
1688 }
1689
1690 // Backward elimination: solves the upper system (U matrix).
1691 for (int i = A.rows() - 1; i >= 0; --i) {
1692
1693 Type sum = Type(0.0);
1694
1695 for (unsigned int j = i + 1; j < A.rows(); ++j)
1696 sum += A(i, j) * b[j];
1697
1698 b[i] = (b[i] - sum) / A(i, i);
1699 }
1700
1701 return b;
1702 }
1703
1704
1712 template<typename Matrix, typename Vector>
1713 inline Vector solve_lu(Matrix A, Vector b) {
1714
1715 // Apply in-place LU decomposition
1717
1718 // Apply forward and backward substitution
1719 return solve_lu_inplace(A, b);
1720 }
1721
1722
1733 template<typename Matrix1, typename Matrix2, typename Vector>
1734 inline Vector solve_lu(const Matrix1& L, const Matrix2& U, const Vector& b) {
1735
1736 Vector x = b;
1737
1738 if (!is_square(L)) {
1739 TH_MATH_ERROR("algebra::solve_lu", L.rows(), MathError::InvalidArgument);
1740 return vec_error(x);
1741 }
1742
1743 if (!is_square(U)) {
1744 TH_MATH_ERROR("algebra::solve_lu", U.rows(), MathError::InvalidArgument);
1745 return vec_error(x);
1746 }
1747
1748 if (L.rows() != U.rows()) {
1749 TH_MATH_ERROR("algebra::solve_lu", U.rows(), MathError::InvalidArgument);
1750 return vec_error(x);
1751 }
1752
1753 if (b.size() != L.rows()) {
1754 TH_MATH_ERROR("algebra::solve_lu", b.size(), MathError::InvalidArgument);
1755 return vec_error(x);
1756 }
1757
1758 using Type1 = matrix_element_t<Matrix1>;
1759
1760 // Forward elimination for L
1761 for (unsigned int i = 1; i < L.rows(); ++i) {
1762
1763 Type1 sum = Type1(0.0);
1764
1765 for (unsigned int j = 0; j < i; ++j)
1766 sum += L(i, j) * x[j];
1767
1768 x[i] -= sum;
1769 }
1770
1771 using Type2 = matrix_element_t<Matrix2>;
1772
1773 // Backward elimination for U
1774 for (int i = U.rows() - 1; i >= 0; --i) {
1775
1776 Type2 sum = Type2(0.0);
1777
1778 for (unsigned int j = i + 1; j < U.rows(); ++j)
1779 sum += U(i, j) * x[j];
1780
1781 if (abs(U(i, i)) < MACH_EPSILON) {
1782 TH_MATH_ERROR("algebra::solve_lu", U(i, i), MathError::DivByZero);
1783 return vec_error(x);
1784 }
1785
1786 x[i] = (x[i] - sum) / U(i, i);
1787 }
1788
1789 return x;
1790 }
1791
1792
1804 template<typename Matrix, typename Vector>
1805 inline Vector solve_cholesky(const Matrix& L, const Vector& b) {
1806
1807 Vector x = b;
1808
1809 if (!is_square(L)) {
1810 TH_MATH_ERROR("algebra::solve_cholesky", L.rows(), MathError::InvalidArgument);
1811 return vec_error(x);
1812 }
1813
1814 if (L.rows() != b.size()) {
1815 TH_MATH_ERROR("algebra::solve_cholesky", b.size(), MathError::InvalidArgument);
1816 return vec_error(x);
1817 }
1818
1819 using Type = matrix_element_t<Matrix>;
1820
1821 // Forward elimination for L
1822 for (unsigned int i = 0; i < L.rows(); ++i) {
1823
1824 Type sum = Type(0.0);
1825
1826 for (unsigned int j = 0; j < i; ++j)
1827 sum += L(i, j) * x[j];
1828
1829 x[i] = (x[i] - sum) / L(i, i);
1830 }
1831
1832 // Backward elimination for L transpose
1833 for (int i = L.rows() - 1; i >= 0; --i) {
1834
1835 Type sum = Type(0.0);
1836
1837 for (unsigned int j = i + 1; j < L.rows(); ++j)
1838 sum += L(j, i) * x[j];
1839
1840 x[i] = (x[i] - sum) / L(i, i);
1841 }
1842
1843 return x;
1844 }
1845
1846
1853 template<typename Matrix, typename Vector>
1854 inline Vector solve(const Matrix& A, const Vector& b) {
1855
1856 // Use LU decomposition to solve the linear system
1857 return solve_lu(A, b);
1858 }
1859
1860
1861 // Other composite operations
1862
1863
1869 template<typename Matrix1, typename Matrix2>
1870 inline Matrix1& inverse(Matrix1& dest, const Matrix2& src) {
1871
1872 if(src.rows() != src.cols()) {
1873 TH_MATH_ERROR("algebra::inverse", src.rows(), MathError::InvalidArgument);
1874 return mat_error(dest);
1875 }
1876
1877 if(dest.rows() != src.rows()) {
1878 TH_MATH_ERROR("algebra::inverse", dest.rows(), MathError::InvalidArgument);
1879 return mat_error(dest);
1880 }
1881
1882 if(dest.cols() != src.cols()) {
1883 TH_MATH_ERROR("algebra::inverse", dest.cols(), MathError::InvalidArgument);
1884 return mat_error(dest);
1885 }
1886
1887 using Type = matrix_element_t<Matrix2>;
1888
1889 // Prepare extended matrix (A|B)
1890 Matrix1 A;
1891 A.resize(src.rows(), src.cols());
1892 dest.resize(src.rows(), src.cols());
1893 mat_copy(A, src);
1894 make_identity(dest);
1895
1896 // Iterate on all columns
1897 for (unsigned int i = 0; i < src.rows(); ++i) {
1898
1899 // Make sure the element on the diagonal
1900 // is non-zero by adding the first non-zero row
1901 if(A(i, i) == (Type) 0) {
1902
1903 bool flag = false;
1904
1905 // Iterate on all rows
1906 for (unsigned int j = i + 1; j < src.rows(); ++j) {
1907
1908 // Add the j-th row to the i-th row
1909 // if Aji is non-zero
1910 if(A(j, i) != (Type) 0) {
1911
1912 for (unsigned int k = 0; k < src.rows(); ++k) {
1913 A(i, k) += A(j, k);
1914 dest(i, k) += dest(j, k);
1915 }
1916
1917 flag = true;
1918 break;
1919 }
1920 }
1921
1922 if(!flag) {
1923 TH_MATH_ERROR("algebra::inverse", flag, MathError::ImpossibleOperation);
1924 return mat_error(dest);
1925 }
1926 }
1927
1928 auto inv_pivot = ((Type) 1.0) / A(i, i);
1929
1930 // Use the current row to make all other
1931 // elements of the column equal to zero
1932 for (unsigned int j = 0; j < src.rows(); ++j) {
1933
1934 // Skip the current row
1935 if(j == i)
1936 continue;
1937
1938 // Multiplication coefficient for
1939 // the elision of Ajk
1940 const auto coeff = A(j, i) * inv_pivot;
1941
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);
1945 }
1946 }
1947
1948 // Divide the current row by the pivot
1949 for (unsigned int j = 0; j < src.rows(); ++j) {
1950 A(i, j) *= inv_pivot;
1951 dest(i, j) *= inv_pivot;
1952 }
1953
1954 }
1955
1956 return dest;
1957 }
1958
1959
1965 template<typename Matrix, typename MatrixInv = Matrix>
1966 inline MatrixInv inverse(const Matrix& m) {
1967 MatrixInv res;
1968 res.resize(m.rows(), m.cols());
1969 inverse(res, m);
1970 return res;
1971 }
1972
1973
1979 template<typename Matrix>
1980 inline Matrix& invert(Matrix& m) {
1981
1982 if(m.rows() != m.cols()) {
1983 TH_MATH_ERROR("algebra::invert", m.rows(), MathError::InvalidArgument);
1984 return mat_error(m);
1985 }
1986
1987 // Prepare extended matrix (A|B)
1988 Matrix tmp;
1989 tmp.resize(m.rows(), m.cols());
1990 inverse(tmp, m);
1991
1992 // Modify the matrix only when the inversion
1993 // has succeeded
1994 mat_copy(m, tmp);
1995 return m;
1996 }
1997
1998
2005 template<typename Matrix>
2006 inline auto det(const Matrix& A) {
2007
2008 Matrix LU = A;
2010
2011 // The determinant of a triangular matrix
2012 // is the product of the elements on its diagonal
2013 return diagonal_product(LU);
2014 }
2015
2016
2017 // Eigensolvers
2018
2019
2027 template<typename Matrix, typename Vector>
2028 inline auto rayleigh_quotient(const Matrix& A, const Vector& x) {
2029
2030 const auto p = dot(x, x);
2031
2032 // Check for division by zero
2033 if (abs(p) < MACH_EPSILON) {
2034 TH_MATH_ERROR("rayleigh_quotient", abs(p), MathError::DivByZero);
2035 return vector_element_t<Vector>(nan());
2036 }
2037
2038 return dot(x, transform(A, x)) / p;
2039 }
2040
2041
2053 template<typename Matrix, typename Vector>
2054 inline auto eigenvalue_power(
2055 const Matrix& A, const Vector& x,
2056 real tolerance = ALGEBRA_EIGEN_TOL, unsigned int max_iter = ALGEBRA_EIGEN_ITER) {
2057
2058 using Type = matrix_element_t<Matrix>;
2059
2060 if (!is_square(A)) {
2061 TH_MATH_ERROR("eigenvalue_power", is_square(A), MathError::InvalidArgument);
2062 return Type(nan());
2063 }
2064
2065 if (x.size() != A.rows()) {
2066 TH_MATH_ERROR("eigenvalue_power", x.size(), MathError::InvalidArgument);
2067 return Type(nan());
2068 }
2069
2070 // Apply a first iteration to initialize
2071 // the current and previous vectors
2072 Vector x_prev = x;
2073 Vector x_curr = normalize(A * x_prev);
2074
2075 // Iteration counter
2076 unsigned int i;
2077
2078 // Iteratively apply the matrix to the vector
2079 for (i = 1; i <= max_iter; ++i) {
2080
2081 x_prev = x_curr;
2082 x_curr = normalize(transform(A, x_prev));
2083
2084 // Stop the algorithm when |x_k+1 +- x_k| is
2085 // less then the tolerance in module
2086 if (norm(x_curr - x_prev) <= tolerance ||
2087 norm(x_curr + x_prev) <= tolerance)
2088 break;
2089 }
2090
2091 // The algorithm did not converge
2092 if (i > max_iter) {
2093 TH_MATH_ERROR("eigenvalue_power", i, MathError::NoConvergence);
2094 return Type(nan());
2095 }
2096
2097 return dot(x_curr, A * x_curr);
2098 }
2099
2100
2114 template<typename Matrix, typename Vector1, typename Vector2 = Vector1>
2115 inline auto eigenpair_power(
2116 const Matrix& A, const Vector1& x, Vector2& v,
2117 real tolerance = ALGEBRA_EIGEN_TOL, unsigned int max_iter = ALGEBRA_EIGEN_ITER) {
2118
2119 using Type = matrix_element_t<Matrix>;
2120
2121 if (!is_square(A)) {
2122 TH_MATH_ERROR("eigenpair_power", is_square(A), MathError::InvalidArgument);
2123 return Type(nan());
2124 }
2125
2126 if (x.size() != A.rows()) {
2127 TH_MATH_ERROR("eigenpair_power", x.size(), MathError::InvalidArgument);
2128 return Type(nan());
2129 }
2130
2131 if (v.size() != x.size()) {
2132 TH_MATH_ERROR("eigenpair_power", v.size(), MathError::InvalidArgument);
2133 return Type(nan());
2134 }
2135
2136 // Apply a first iteration to initialize
2137 // the current and previous vectors
2138 Vector1 x_prev = x;
2139 Vector1 x_curr = normalize(A * x_prev);
2140
2141 // Iteration counter
2142 unsigned int i;
2143
2144 // Iteratively apply the matrix to the vector
2145 for (i = 1; i <= max_iter; ++i) {
2146
2147 x_prev = x_curr;
2148 x_curr = normalize(transform(A, x_prev));
2149
2150 // Stop the algorithm when |x_k+1 +- x_k| is
2151 // less then the tolerance in module
2152 if (norm(x_curr - x_prev) <= tolerance ||
2153 norm(x_curr + x_prev) <= tolerance)
2154 break;
2155 }
2156
2157 // The algorithm did not converge
2158 if (i > max_iter) {
2159 TH_MATH_ERROR("eigenpair_power", i, MathError::NoConvergence);
2160 return Type(nan());
2161 }
2162
2163 // Overwrite with the eigenvector
2164 vec_copy(v, x_curr);
2165
2166 return dot(x_curr, A * x_curr);
2167 }
2168
2169
2181 template<typename Matrix, typename Vector>
2183 const Matrix& A, const Vector& x,
2184 real tolerance = ALGEBRA_EIGEN_TOL,
2185 unsigned int max_iter = ALGEBRA_EIGEN_ITER) {
2186
2187 using Type = matrix_element_t<Matrix>;
2188
2189 if (!is_square(A)) {
2190 TH_MATH_ERROR("eigenvalue_inverse", is_square(A), MathError::InvalidArgument);
2191 return Type(nan());
2192 }
2193
2194 if (A.rows() != x.size()) {
2195 TH_MATH_ERROR("eigenvalue_inverse", x.size(), MathError::InvalidArgument);
2196 return Type(nan());
2197 }
2198
2199 // Compute the LU decomposition of A to speed up system solution
2200 Matrix LU = A;
2202
2203 // Compute the first step to initialize the two vectors
2204 Vector x_prev = normalize(x);
2205 Vector x_curr = x_prev;
2206 solve_lu_inplace(LU, x_curr);
2207
2208 // Iteration counter
2209 unsigned int i;
2210
2211 for (i = 1; i <= max_iter; ++i) {
2212
2213 // Solve the linear system using the LU decomposition
2214 x_prev = x_curr;
2215 make_normalized(solve_lu_inplace(LU, x_curr));
2216
2217 // Stop the algorithm when |x_k+1 +- x_k| is
2218 // less then the tolerance in module
2219 if (norm(x_curr - x_prev) <= tolerance ||
2220 norm(x_curr + x_prev) <= tolerance)
2221 break;
2222 }
2223
2224 // The algorithm did not converge
2225 if (i > max_iter) {
2226 TH_MATH_ERROR("eigenvalue_inverse", i, MathError::NoConvergence);
2227 return Type(nan());
2228 }
2229
2230 // A and inverse(A) have the same eigenvectors
2231 return dot(x_curr, A * x_curr);
2232 }
2233
2234
2248 template<typename Matrix, typename Vector1, typename Vector2>
2250 const Matrix& A, const Vector1& x, Vector2& v,
2251 real tolerance = ALGEBRA_EIGEN_TOL, unsigned int max_iter = ALGEBRA_EIGEN_ITER) {
2252
2253 using Type = matrix_element_t<Matrix>;
2254
2255 if (!is_square(A)) {
2256 TH_MATH_ERROR("eigenpair_inverse", is_square(A), MathError::InvalidArgument);
2257 return Type(nan());
2258 }
2259
2260 if (A.rows() != x.size()) {
2261 TH_MATH_ERROR("eigenpair_inverse", x.size(), MathError::InvalidArgument);
2262 return Type(nan());
2263 }
2264
2265 // Compute the LU decomposition of A to speed up system solution
2266 Matrix LU = A;
2268
2269 // Compute the first step to initialize the two vectors
2270 Vector1 x_prev = normalize(x);
2271 Vector1 x_curr = x_prev;
2272 solve_lu_inplace(LU, x_curr);
2273
2274 // Iteration counter
2275 unsigned int i;
2276
2277 for (i = 1; i <= max_iter; ++i) {
2278
2279 // Solve the linear system using the LU decomposition
2280 x_prev = x_curr;
2281 make_normalized(solve_lu_inplace(LU, x_curr));
2282
2283 // Stop the algorithm when |x_k+1 +- x_k| is
2284 // less then the tolerance in module
2285 if (norm(x_curr - x_prev) <= tolerance ||
2286 norm(x_curr + x_prev) <= tolerance)
2287 break;
2288 }
2289
2290 // The algorithm did not converge
2291 if (i > max_iter) {
2292 TH_MATH_ERROR("eigenpair_inverse", i, MathError::NoConvergence);
2293 return Type(nan());
2294 }
2295
2296 // Overwrite with the eigenvector
2297 vec_copy(v, x_curr);
2298
2299 // A and inverse(A) have the same eigenvectors
2300 return dot(x_curr, A * x_curr);
2301 }
2302
2303
2320 template<typename Matrix, typename Vector, typename T = matrix_element_t<Matrix>>
2321 inline Vector eigenvector_inverse(
2322 const Matrix& A, const T& lambda, const Vector& x,
2323 real tolerance = ALGEBRA_EIGEN_TOL, unsigned int max_iter = ALGEBRA_EIGEN_ITER) {
2324
2325 if (!is_square(A)) {
2326 TH_MATH_ERROR("eigenvector_inverse", is_square(A), MathError::InvalidArgument);
2327 return make_error<Vector>(A.rows());
2328 }
2329
2330 if (A.rows() != x.size()) {
2331 TH_MATH_ERROR("eigenvector_inverse", x.size(), MathError::InvalidArgument);
2332 return make_error<Vector>(A.rows());
2333 }
2334
2335 Matrix LU = A;
2336 mat_shift_diagonal(LU, -lambda);
2337
2338 // Compute the LU decomposition of A
2339 // to speed up system solution
2341
2342 // Compute the first step to initialize the two vectors
2343 Vector v_prev = normalize(x);
2344 Vector v_curr = v_prev;
2345 solve_lu_inplace(LU, v_curr);
2346
2347 // Iteration counter
2348 unsigned int i;
2349
2350 for (i = 1; i <= max_iter; ++i) {
2351
2352 // Solve the linear system using the LU decomposition
2353 v_prev = v_curr;
2354 make_normalized(solve_lu_inplace(LU, v_curr));
2355
2356 // Stop the algorithm when |x_k+1 +- x_k| is
2357 // less then the tolerance in module
2358 if (norm(v_curr - v_prev) <= tolerance ||
2359 norm(v_curr + v_prev) <= tolerance)
2360 break;
2361 }
2362
2363 // The algorithm did not converge
2364 if (i > max_iter) {
2365 TH_MATH_ERROR("eigenvector_inverse", i, MathError::NoConvergence);
2366 return vec_error(v_curr);
2367 }
2368
2369 return v_curr;
2370 }
2371
2372
2386 template<typename Matrix, typename Vector, typename T = matrix_element_t<Matrix>>
2388 const Matrix& A, const T& lambda, const Vector& x,
2389 real tolerance = ALGEBRA_EIGEN_TOL, unsigned int max_iter = ALGEBRA_EIGEN_ITER) {
2390
2391 using Type = matrix_element_t<Matrix>;
2392
2393 if (!is_square(A)) {
2394 TH_MATH_ERROR("eigenvalue_rayleigh", is_square(A), MathError::InvalidArgument);
2395 return Type(nan());
2396 }
2397
2398 if (A.rows() != x.size()) {
2399 TH_MATH_ERROR("eigenvalue_rayleigh", x.size(), MathError::InvalidArgument);
2400 return Type(nan());
2401 }
2402
2403 // Keep track of the shifted matrix
2404 Matrix A_shift = A;
2405 mat_shift_diagonal(A_shift, -lambda);
2406
2407 Type lambda_prev = lambda;
2408 Type lambda_curr = lambda;
2409
2410 Vector x_prev = normalize(x);
2411 Vector x_curr = normalize(solve(A_shift, x));
2412
2413 unsigned int i;
2414
2415 for (i = 1; i <= max_iter; ++i) {
2416
2417 // Update the eigenvalue approximation
2418 // using Rayleigh quotient (avoiding normalization)
2419 lambda_prev = lambda_curr;
2420 lambda_curr = dot(x_curr, transform(A, x_curr));
2421
2422 // Shift the diagonal by the difference between
2423 // subsequent eigenvalues steps, to avoid copying matrix A
2424 mat_shift_diagonal(A_shift, lambda_prev - lambda_curr);
2425
2426 // Solve the linear system each time, as the shifted
2427 // matrix is continuously updated
2428 x_prev = x_curr;
2429 x_curr = normalize(solve(A_shift, x_curr));
2430
2431 // Stop the algorithm when |x_k+1 +- x_k| is
2432 // less then the tolerance in module
2433 if (norm(x_curr - x_prev) <= tolerance ||
2434 norm(x_curr + x_prev) <= tolerance)
2435 break;
2436 }
2437
2438 if (i > max_iter) {
2439 TH_MATH_ERROR("eigenvalue_rayleigh", i, MathError::NoConvergence);
2440 return Type(nan());
2441 }
2442
2443 return dot(x_curr, transform(A, x_curr));
2444 }
2445
2446
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,
2466 real tolerance = ALGEBRA_EIGEN_TOL, unsigned int max_iter = ALGEBRA_EIGEN_ITER) {
2467
2468 using Type = matrix_element_t<Matrix>;
2469
2470 if (!is_square(A)) {
2471 TH_MATH_ERROR("eigenpair_rayleigh", is_square(A), MathError::InvalidArgument);
2472 return Type(nan());
2473 }
2474
2475 if (A.rows() != x.size()) {
2476 TH_MATH_ERROR("eigenpair_rayleigh", x.size(), MathError::InvalidArgument);
2477 return Type(nan());
2478 }
2479
2480 // Keep track of the shifted matrix
2481 Matrix A_shift = A;
2482 mat_shift_diagonal(A_shift, -lambda);
2483
2484 Type lambda_prev = lambda;
2485 Type lambda_curr = lambda;
2486
2487 Vector1 x_prev = normalize(x);
2488 Vector1 x_curr = normalize(solve(A_shift, x));
2489
2490 unsigned int i;
2491
2492 for (i = 1; i <= max_iter; ++i) {
2493
2494 // Update the eigenvalue approximation
2495 // using Rayleigh quotient (avoiding normalization)
2496 lambda_prev = lambda_curr;
2497 lambda_curr = dot(x_curr, transform(A, x_curr));
2498
2499 // Shift the diagonal by the difference between
2500 // subsequent eigenvalues steps, to avoid copying matrix A
2501 mat_shift_diagonal(A_shift, lambda_prev - lambda_curr);
2502
2503 // Solve the linear system each time, as the shifted
2504 // matrix is continuously updated
2505 x_prev = x_curr;
2506 x_curr = normalize(solve(A_shift, x_curr));
2507
2508 // Stop the algorithm when |x_k+1 +- x_k| is
2509 // less then the tolerance in module
2510 if (norm(x_curr - x_prev) <= tolerance ||
2511 norm(x_curr + x_prev) <= tolerance)
2512 break;
2513 }
2514
2515 if (i > max_iter) {
2516 TH_MATH_ERROR("eigenpair_rayleigh", i, MathError::NoConvergence);
2517 return Type(nan());
2518 }
2519
2520 // Overwrite with the eigenvector
2521 vec_copy(v, x_curr);
2522
2523 return dot(x_curr, transform(A, x_curr));
2524 }
2525 }
2526}
2527
2528#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: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