Easy3D 2.5.3
mat.h
1/********************************************************************
2 * Copyright (C) 2015 Liangliang Nan <liangliang.nan@gmail.com>
3 * https://3d.bk.tudelft.nl/liangliang/
4 *
5 * This file is part of Easy3D. If it is useful in your research/work,
6 * I would be grateful if you show your appreciation by citing it:
7 * ------------------------------------------------------------------
8 * Liangliang Nan.
9 * Easy3D: a lightweight, easy-to-use, and efficient C++ library
10 * for processing and rendering 3D data.
11 * Journal of Open Source Software, 6(64), 3255, 2021.
12 * ------------------------------------------------------------------
13 *
14 * Easy3D is free software; you can redistribute it and/or modify
15 * it under the terms of the GNU General Public License Version 3
16 * as published by the Free Software Foundation.
17 *
18 * Easy3D is distributed in the hope that it will be useful,
19 * but WITHOUT ANY WARRANTY; without even the implied warranty of
20 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 * GNU General Public License for more details.
22 *
23 * You should have received a copy of the GNU General Public License
24 * along with this program. If not, see <http://www.gnu.org/licenses/>.
25 ********************************************************************/
26
27#ifndef EASY3D_CORE_MAT_H
28#define EASY3D_CORE_MAT_H
29
30#include <iostream>
31#include <iomanip>
32
33#include <easy3d/core/vec.h>
34#include <easy3d/core/constant.h>
35#include <easy3d/util/logging.h>
36
37
38namespace easy3d {
39
40
44 //#define MATRIX_ROW_MAJOR
45
46 template <typename T> class Mat2;
47 template <typename T> class Mat3;
48 template <typename T> class Mat4;
49
50
64 template <size_t N, size_t M, typename T>
65 class Mat
66 {
67 public:
68 // ----------------- constructor and destructor -------------------
69
75 Mat() = default;
76
78 explicit Mat(T s);
79
92 template <size_t rN, size_t rM>
93 explicit Mat(const Mat<rN, rM, T> &rhs);
94
96 explicit Mat(const T *m);
97
101
102 // ------------------- size query ---------------------
103
105 size_t num_rows() const { return N; }
106
108 size_t num_columns() const { return M; }
109
110 // --------------------- access -----------------------
111
113 Vec<M, T> row(size_t r) const;
114
116 Vec<N, T> col(size_t c) const;
117
119 const T& operator()(size_t r, size_t c) const;
120
122 T& operator()(size_t r, size_t c);
123
125 operator const T*() const;
126
128 operator T*();
129
130 // --------------------- modification ------------------
131
133 void load_zero();
134
136 void load_identity(T s = T(1));
137
143 template <size_t vN>
144 void set_row(size_t r, const Vec<vN, T> &v);
145
151 template <size_t vN>
152 void set_col(size_t c, const Vec<vN, T> &v);
153
155 void swap_rows(size_t a, size_t b);
156
158 void swap_cols(size_t a, size_t b);
159
160 // ------- matrix-matrix arithmetic operators --------
161
163 bool operator==(const Mat<N, M, T> &rhs) const;
164
166 bool operator!=(const Mat<N, M, T> &rhs) const;
167
173 template <size_t rM>
175
178
181
184
185 // ------- matrix-vector arithmetic operators --------
186
191 Vec<N, T> operator*(const Vec<M, T> &rhs) const;
192
193 // ------- matrix-scalar arithmetic operators --------
194
199
200 // ------- matrix-matrix assignment operators --------
201
204
207
210
211 // ------- matrix-scalar assignment operators --------
212
215
218
221
224
225 protected:
227 T m_[N * M];
228 };
229
230
231 // ------- global scalar-matrix arithmetic operators --------
232
234 template <size_t N, size_t M, typename T>
235 Mat<N, M, T> operator*(T lhs, const Mat<N, M, T> &rhs);
236
237 // ------- global matrix-matrix multiplication operators --------
238
240 template <typename T>
241 Mat2<T> operator*(const Mat2<T> &lhs, const Mat2<T> &rhs);
242
244 template <typename T>
245 Mat3<T> operator*(const Mat3<T> &lhs, const Mat3<T> &rhs);
246
248 template <typename T>
249 Mat4<T> operator*(const Mat4<T> &lhs, const Mat4<T> &rhs);
250
251 // ------- global matrix-vector arithmetic operators --------
252
257 template <typename T>
258 Vec<3, T> operator*(const Mat4<T> &lhs, const Vec<3, T> &rhs);
259
264 template <typename T>
265 Vec<2, T> operator*(const Mat3<T> &lhs, const Vec<2, T> &rhs);
266
270 template <typename T>
271 Vec<2, T> operator*(const Mat2<T> &lhs, const Vec<2, T> &rhs);
275 template <typename T>
276 Vec<3, T> operator*(const Mat3<T> &lhs, const Vec<3, T> &rhs);
280 template <typename T>
281 Vec<4, T> operator*(const Mat4<T> &lhs, const Vec<4, T> &rhs);
282
283 // -------------- global matrix related function ---------------
284
286 template <size_t N, typename T>
287 T trace(const Mat<N, N, T> &m);
288
294 template <size_t N, typename T, size_t A>
296
298 template <size_t N, size_t M, typename T>
300
306 template <size_t N, typename T>
308
313 template <size_t M, size_t N, typename T>
314 Mat<N, M, T> tensor(const Vec<M, T> &u, const Vec<N, T> &v);
315
325 template<size_t N, size_t M, typename T>
327
341 template<size_t N, typename T>
342 bool lu_decomposition(const Mat<N, N, T> &a, Mat<N, N, T> *alu, Vec<N, T> *rowp, T *d);
343
401 template <size_t N, typename T>
402 void lu_back_substitution(const Mat<N, N, T> &alu, const Vec<N, T> &rowp, const Vec<N, T> &b, Vec<N, T> *x);
403
411 template<size_t N, typename FT>
413
420 template <size_t N, typename FT>
421 void cholesky_solve(const Mat<N, N, FT> &L, const Vec<N, FT> &b, Vec<N, FT>& x);
422
430 template<size_t N, size_t M, typename T>
431 void cholesky_solve(const Mat<N, N, T> &L, const Mat<N, M, T> &B, Mat<N, M, T>& X);
432
433
434 template <size_t N, size_t M, typename T>
435 std::ostream& operator<< (std::ostream& output, const Mat<N, M, T>& m);
436
437 template <size_t N, size_t M, typename T>
438 std::istream& operator>> (std::istream& input, Mat<N, M, T>& m);
439
440
441
442 /*******************************************************************************
443
444 IMPLEMENTATION:
445
446 *******************************************************************************/
447
448 /*----------------------------------------------------------------------------*/
449 template <size_t N, size_t M, typename T>
450 inline Mat<N, M, T>::Mat(T s) {
451 for (size_t i = 0; i < N; ++i) {
452 for (size_t j = 0; j < M; ++j)
453 if (i == j)
454 (*this)(i, j) = T(s);
455 else
456 (*this)(i, j) = T(0);
457 }
458 }
459
460 /*----------------------------------------------------------------------------*/
461 template <size_t N, size_t M, typename T>
462 template <size_t rN, size_t rM>
464 assert(rN >= N);
465 assert(rM >= M);
466
467 for (size_t i = 0; i < N; ++i)
468 for (size_t j = 0; j < M; ++j)
469 (*this)(i, j) = rhs(i, j);
470 }
471
472 /*----------------------------------------------------------------------------*/
473 template <size_t N, size_t M, typename T>
474 inline Mat<N, M, T>::Mat(const T *m) {
475 assert(m != 0);
476
477 for (size_t i = 0; i < N * M; ++i)
478 m_[i] = m[i];
479 }
480
481 /*----------------------------------------------------------------------------*/
482 template <size_t N, size_t M, typename T>
484 Mat<N, M, T> result;
485 for (size_t i = 0; i < N; ++i) {
486 for (size_t j = 0; j < M; ++j)
487 if (i == j)
488 result(i, j) = T(1);
489 else
490 result(i, j) = T(0);
491 }
492 return result;
493 }
494
495
496 /*----------------------------------------------------------------------------*/
497 template <size_t N, size_t M, typename T>
498 inline Vec<M, T> Mat<N, M, T>::row(size_t r) const {
499 assert(r < N);
500
501 Vec<M, T> result;
502 for (size_t i = 0; i < M; ++i)
503 result[i] = (*this)(r, i);
504 return result;
505 }
506
507 /*----------------------------------------------------------------------------*/
508 template <size_t N, size_t M, typename T>
509 inline Vec<N, T> Mat<N, M, T>::col(size_t c) const {
510 assert(c < M);
511
512 Vec<N, T> result;
513 for (size_t i = 0; i < N; ++i)
514 result[i] = (*this)(i, c);
515 return result;
516 }
517
518 /*----------------------------------------------------------------------------*/
519 template <size_t N, size_t M, typename T>
520 template <size_t vN>
521 inline void Mat<N, M, T>::set_row(size_t r, const Vec<vN, T> &v) {
522 assert(r < N);
523 assert(vN >= M);
524
525 for (size_t i = 0; i < M; ++i)
526 (*this)(r, i) = v[i];
527 }
528
529 /*----------------------------------------------------------------------------*/
530 template <size_t N, size_t M, typename T>
531 template <size_t vN>
532 inline void Mat<N, M, T>::set_col(size_t c, const Vec<vN, T> &v) {
533 assert(c < M);
534 assert(vN >= N);
535
536 for (size_t i = 0; i < N; ++i)
537 (*this)(i, c) = v[i];
538 }
539
540 /*----------------------------------------------------------------------------*/
541 template <size_t N, size_t M, typename T>
542 inline void Mat<N, M, T>::swap_rows(size_t a, size_t b) {
543 assert(a < N);
544 assert(b < N);
545
546 for (size_t i = 0; i < M; ++i) {
547 T tmp = (*this)(a, i);
548 (*this)(a, i) = (*this)(b, i);
549 (*this)(b, i) = tmp;
550 }
551 }
552
553 /*----------------------------------------------------------------------------*/
554 template <size_t N, size_t M, typename T>
556 size_t size = N * M;
557 for (size_t i = 0; i < size; ++i)
558 m_[i] = T(0);
559 }
560
561 /*----------------------------------------------------------------------------*/
562 template <size_t N, size_t M, typename T>
563 inline void Mat<N, M, T>::load_identity(T s /* = T(1)*/) {
564 for (size_t i = 0; i < N; ++i) {
565 for (size_t j = 0; j < M; ++j)
566 if (i == j)
567 (*this)(i, j) = T(s);
568 else
569 (*this)(i, j) = T(0);
570 }
571 }
572
573 /*----------------------------------------------------------------------------*/
574 template <size_t N, size_t M, typename T>
575 inline void Mat<N, M, T>::swap_cols(size_t a, size_t b) {
576 assert(a < M);
577 assert(b < M);
578
579 for (size_t i = 0; i < N; ++i) {
580 T tmp = (*this)(i, a);
581 (*this)(i, a) = (*this)(i, b);
582 (*this)(i, b) = tmp;
583 }
584 }
585
586 /*----------------------------------------------------------------------------*/
587 template <size_t N, size_t M, typename T>
588 inline Mat<N, M, T>::operator const T*() const {
589 return m_;
590 }
591
592 /*----------------------------------------------------------------------------*/
593 template <size_t N, size_t M, typename T>
595 return m_;
596 }
597
598 /*----------------------------------------------------------------------------*/
599 template <size_t N, size_t M, typename T>
600 inline const T& Mat<N, M, T>::operator()(size_t row, size_t col) const {
601 assert(row < N);
602 assert(col < M);
603
604 #ifdef MATRIX_ROW_MAJOR
605 return m_[row * M + col]; // row-major
606 #else
607 return m_[col * N + row]; // column-major
608 #endif
609 }
610
611 /*----------------------------------------------------------------------------*/
612 template <size_t N, size_t M, typename T>
613 inline T& Mat<N, M, T>::operator()(size_t row, size_t col) {
614 assert(row < N);
615 assert(col < M);
616
617 #ifdef MATRIX_ROW_MAJOR
618 return m_[row * M + col]; // row-major
619 #else
620 return m_[col * N + row]; // column-major
621 #endif
622 }
623
624 /*----------------------------------------------------------------------------*/
625 template <size_t N, size_t M, typename T>
626 inline bool Mat<N, M, T>::operator==(const Mat<N, M, T> &rhs) const {
627 bool result = true;
628 for (size_t i = 0; i < N * M; ++i)
629 result &= equal(m_[i], rhs[i]);
630 return result;
631 }
632
633 /*----------------------------------------------------------------------------*/
634 template <size_t N, size_t M, typename T>
635 inline bool Mat<N, M, T>::operator!=(const Mat<N, M, T> &rhs) const {
636 return !(*this == rhs);
637 }
638
639
640 /*----------------------------------------------------------------------------*/
641 template <size_t N, size_t M, typename T>
642 template <size_t rM>
644 Mat<N, rM, T> result;
645 for (size_t i = 0; i < N; ++i) {
646 for (size_t j = 0; j < rM; ++j) {
647 result(i, j) = T(0);
648 for (size_t k = 0; k < M; ++k) {
649 result(i, j) += (*this)(i, k) * rhs(k, j);
650 }
651 }
652 }
653 return result;
654 }
655
656 /*----------------------------------------------------------------------------*/
657 template <size_t N, size_t M, typename T>
659 Mat<N, M, T> result;
660 for (size_t i = 0; i < N * M; ++i)
661 result[i] = m_[i] + rhs[i];
662 return result;
663 }
664
665 /*----------------------------------------------------------------------------*/
666 template <size_t N, size_t M, typename T>
668 Mat<N, M, T> result;
669 for (size_t i = 0; i < N * M; ++i)
670 result[i] = m_[i] - rhs[i];
671 return result;
672 }
673
674 /*----------------------------------------------------------------------------*/
675 template <size_t N, size_t M, typename T>
677 Mat<N, M, T> result;
678 for (size_t i = 0; i < N * M; ++i)
679 result[i] = -m_[i];
680 return result;
681 }
682
683
684 // MATRIX-VECTOR ARITHMETIC OPERATORS:
685
686 /*----------------------------------------------------------------------------*/
687 template <size_t N, size_t M, typename T>
688 inline Vec<N, T> Mat<N, M, T>::operator*(const Vec<M, T> &rhs) const {
689 Vec<N, T> result;
690 for (size_t i = 0; i < N; ++i) {
691 result[i] = 0;
692 for (size_t j = 0; j < M; ++j) {
693 result[i] += (*this)(i, j) * rhs[j];
694 }
695 }
696 return result;
697 }
698
699
700 // MATRIX-SCALAR ARITHMETIC OPERATORS:
701
702 /*----------------------------------------------------------------------------*/
703 template <size_t N, size_t M, typename T>
705 Mat<N, M, T> result;
706 for (size_t i = 0; i < N * M; ++i)
707 result[i] = m_[i] * rhs;
708 return result;
709 }
710
711 /*----------------------------------------------------------------------------*/
712 template <size_t N, size_t M, typename T>
714 Mat<N, M, T> result;
715 for (size_t i = 0; i < N * M; ++i)
716 result[i] = m_[i] / rhs;
717 return result;
718 }
719
720 // MATRIX-MATRIX ASSIGNMENT OPERATORS:
721
722 /*----------------------------------------------------------------------------*/
723 template <size_t N, size_t M, typename T>
725 Mat<N, M, T> tmp;
726 tmp = *this * rhs;
727 *this = tmp;
728 return *this;
729 }
730
731 /*----------------------------------------------------------------------------*/
732 template <size_t N, size_t M, typename T>
734 for (size_t i = 0; i < N * M; ++i)
735 m_[i] += rhs[i];
736 return *this;
737 }
738
739 /*----------------------------------------------------------------------------*/
740 template <size_t N, size_t M, typename T>
742 for (size_t i = 0; i < N * M; ++i)
743 m_[i] -= rhs[i];
744 return *this;
745 }
746
747
748 // MATRIX-SCALAR ASSIGNMENT OPERATORS:
749
750 /*----------------------------------------------------------------------------*/
751 template <size_t N, size_t M, typename T>
753 for (size_t i = 0; i < N * M; ++i)
754 m_[i] *= rhs;
755 return *this;
756 }
757
758 /*----------------------------------------------------------------------------*/
759 template <size_t N, size_t M, typename T>
761 for (size_t i = 0; i < N * M; ++i)
762 m_[i] /= rhs;
763 return *this;
764 }
765
766 /*----------------------------------------------------------------------------*/
767 template <size_t N, size_t M, typename T>
769 for (size_t i = 0; i < N * M; ++i)
770 m_[i] += rhs;
771 return *this;
772 }
773
774 /*----------------------------------------------------------------------------*/
775 template <size_t N, size_t M, typename T>
777 for (size_t i = 0; i < N * M; ++i)
778 m_[i] -= rhs;
779 return *this;
780 }
781
782
783 // GLOBAL SCALAR-MATRIX ARITHMETIC OPERATORS:
784
785 /*----------------------------------------------------------------------------*/
786 template <size_t N, size_t M, typename T>
787 inline Mat<N, M, T> operator*(T lhs, const Mat<N, M, T> &rhs) {
788 Mat<N, M, T> result;
789 for (size_t i = 0; i < N * M; ++i)
790 result[i] = lhs * rhs[i];
791 return result;
792 }
793
794
795 // ------- global matrix-matrix multiplication operators --------
796
797 template <typename T>
798 inline Mat2<T> operator*(const Mat2<T> &lhs, const Mat2<T> &rhs) {
799 Mat2<T> result;
800 for (size_t i = 0; i < 2; ++i) {
801 for (size_t j = 0; j < 2; ++j) {
802 result(i, j) = T(0);
803 for (size_t k = 0; k < 2; ++k) {
804 result(i, j) += lhs(i, k) * rhs(k, j);
805 }
806 }
807 }
808 return result;
809 }
810
811 template <typename T>
812 inline Mat3<T> operator*(const Mat3<T> &lhs, const Mat3<T> &rhs) {
813 Mat3<T> result;
814 for (size_t i = 0; i < 3; ++i) {
815 for (size_t j = 0; j < 3; ++j) {
816 result(i, j) = T(0);
817 for (size_t k = 0; k < 3; ++k) {
818 result(i, j) += lhs(i, k) * rhs(k, j);
819 }
820 }
821 }
822 return result;
823 }
824
825 template <typename T>
826 inline Mat4<T> operator*(const Mat4<T> &lhs, const Mat4<T> &rhs) {
827 Mat4<T> result;
828 for (size_t i = 0; i < 4; ++i) {
829 for (size_t j = 0; j < 4; ++j) {
830 result(i, j) = T(0);
831 for (size_t k = 0; k < 4; ++k) {
832 result(i, j) += lhs(i, k) * rhs(k, j);
833 }
834 }
835 }
836 return result;
837 }
838
839 // GLOBAL MATRIX-VECTOR ARITHMETIC OPERATORS:
840
841 /*----------------------------- homogeneous version --------------------------*/
842 template <typename T>
843 inline Vec<3, T> operator*(const Mat4<T> &lhs, const Vec<3, T> &rhs) {
844 Vec<4, T> tmp(rhs, T(1));
845 Vec<4, T> result = lhs * tmp;
846 result /= result[3];
847 return Vec<3, T>((T*)result);
848 }
849
850 /*----------------------------- homogeneous version --------------------------*/
851 template <typename T>
852 inline Vec<2, T> operator*(const Mat3<T> &lhs, const Vec<2, T> &rhs) {
853 Vec<3, T> tmp(rhs, T(1));
854 Vec<3, T> result = lhs * tmp;
855 result /= result[2];
856 return Vec<2, T>((T*)result);
857 }
858
859
860 /*--------------------------- non-homogeneous version ------------------------*/
861 template <typename T>
862 inline Vec<2, T> operator*(const Mat2<T> &lhs, const Vec<2, T> &rhs) {
863 Vec<2, T> result;
864 for (size_t i = 0; i < 2; ++i) {
865 result[i] = 0;
866 for (size_t j = 0; j < 2; ++j) {
867 result[i] += lhs(i, j) * rhs[j];
868 }
869 }
870 return result;
871 }
872
873 /*--------------------------- non-homogeneous version ------------------------*/
875 template <typename T>
876 inline Vec<3, T> operator*(const Mat3<T> &lhs, const Vec<3, T> &rhs) {
877 Vec<3, T> result;
878 for (size_t i = 0; i < 3; ++i) {
879 result[i] = 0;
880 for (size_t j = 0; j < 3; ++j) {
881 result[i] += lhs(i, j) * rhs[j];
882 }
883 }
884 return result;
885 }
886
887 /*--------------------------- non-homogeneous version ------------------------*/
889 template <typename T>
890 inline Vec<4, T> operator*(const Mat4<T> &lhs, const Vec<4, T> &rhs) {
891 Vec<4, T> result;
892 for (size_t i = 0; i < 4; ++i) {
893 result[i] = 0;
894 for (size_t j = 0; j < 4; ++j) {
895 result[i] += lhs(i, j) * rhs[j];
896 }
897 }
898 return result;
899 }
900
901
902
903 // GLOBAL SERVICES:
904
905 /*----------------------------------------------------------------------------*/
906 template <size_t N, size_t M, typename T>
908 Mat<M, N, T> result;
909 for (size_t i = 0; i < N; ++i) {
910 for (size_t j = 0; j < M; ++j) {
911 result(j, i) = m(i, j);
912 }
913 }
914 return result;
915 }
916
917 /*----------------------------------------------------------------------------*/
919 template <size_t D, typename T>
920 inline T trace(const Mat<D, D, T> &m) {
921 T result = m(0, 0);
922 for (size_t i = 1; i < D; ++i)
923 result += m(i, i);
924 return result;
925 }
926
927 /*----------------------------------------------------------------------------*/
929 template <size_t N, typename T>
930 inline T determinant(const Mat<N, N, T> &m) {
931 /* Use LU decomposition to find the determinant. */
932 Mat<N, N, T> tmp;
933 Vec<N, size_t> rowp;
934 T result;
935 lu_decomposition(m, &tmp, &rowp, &result);
936 for (size_t i = 0; i < N; ++i)
937 result *= tmp(i, i);
938 return result;
939 }
940
942 template <typename T>
943 inline T determinant(const Mat2<T> &m) {
944 return m(0, 0) * m(1, 1) - m(0, 1) * m(1, 0);
945 }
946
948 template <typename T>
949 inline T determinant(const Mat3<T> &m) {
950 return
951 m(0, 0) * (m(1, 1) * m(2, 2) - m(2, 1) * m(1, 2)) +
952 m(0, 1) * (m(2, 0) * m(1, 2) - m(1, 0) * m(2, 2)) +
953 m(0, 2) * (m(1, 0) * m(2, 1) - m(2, 0) * m(1, 1));
954 }
955
957 template <typename T>
958 inline T determinant(const Mat4<T> &m) {
959 return
960 m(0, 3) * m(1, 2) * m(2, 1) * m(3, 0) - m(0, 2) * m(1, 3) * m(2, 1) * m(3, 0) -
961 m(0, 3) * m(1, 1) * m(2, 2) * m(3, 0) + m(0, 1) * m(1, 3) * m(2, 2) * m(3, 0) +
962 m(0, 2) * m(1, 1) * m(2, 3) * m(3, 0) - m(0, 1) * m(1, 2) * m(2, 3) * m(3, 0) -
963 m(0, 3) * m(1, 2) * m(2, 0) * m(3, 1) + m(0, 2) * m(1, 3) * m(2, 0) * m(3, 1) +
964 m(0, 3) * m(1, 0) * m(2, 2) * m(3, 1) - m(0, 0) * m(1, 3) * m(2, 2) * m(3, 1) -
965 m(0, 2) * m(1, 0) * m(2, 3) * m(3, 1) + m(0, 0) * m(1, 2) * m(2, 3) * m(3, 1) +
966 m(0, 3) * m(1, 1) * m(2, 0) * m(3, 2) - m(0, 1) * m(1, 3) * m(2, 0) * m(3, 2) -
967 m(0, 3) * m(1, 0) * m(2, 1) * m(3, 2) + m(0, 0) * m(1, 3) * m(2, 1) * m(3, 2) +
968 m(0, 1) * m(1, 0) * m(2, 3) * m(3, 2) - m(0, 0) * m(1, 1) * m(2, 3) * m(3, 2) -
969 m(0, 2) * m(1, 1) * m(2, 0) * m(3, 3) + m(0, 1) * m(1, 2) * m(2, 0) * m(3, 3) +
970 m(0, 2) * m(1, 0) * m(2, 1) * m(3, 3) - m(0, 0) * m(1, 2) * m(2, 1) * m(3, 3) -
971 m(0, 1) * m(1, 0) * m(2, 2) * m(3, 3) + m(0, 0) * m(1, 1) * m(2, 2) * m(3, 3);
972 }
973
974 /*----------------------------------------------------------------------------*/
976 template <size_t N, typename T>
978 // Use Gauss-Jordan elimination to find inverse. This uses less memory than the LU decomposition method.
979 // See lu_back_substitution() for an example of computing inverse using LU decomposition.
980 size_t indxc[N], indxr[N], ipiv[N];
981
982 Mat<N, N, T> result = m;
983
984 for (size_t i = 0; i < N; ++i)
985 ipiv[i] = 0;
986
987 for (size_t i = 0; i < N; ++i) { // for each column
988 T max = T(0);
989 size_t maxc, maxr;
990 for (size_t j = 0; j < N; ++j) { // search for pivot
991 if (ipiv[j] != 1) {
992 for (size_t k = 0; k < N; ++k) { // for each column
993 if (ipiv[k] == 0) {
994 T element = std::abs(result(j, k));
995 if (element > max) {
996 max = element;
997 maxr = j;
998 maxc = k;
999 }
1000 }
1001 }
1002 }
1003 }
1004 ++ipiv[maxc];
1005
1006 if (maxr != maxc)
1007 result.swap_rows(maxr, maxc);
1008
1009 indxr[i] = maxr;
1010 indxc[i] = maxc;
1011
1012 // check for singular matrix:
1013 if (std::abs(result(maxc, maxc)) < epsilon<T>()) {
1014 LOG(ERROR) << "input matrix is singular";
1015 return result; // return partial result
1016 }
1017
1018 // multiply row by 1/pivot:
1019 T rpivot = T(1) / result(maxc, maxc);
1020 result(maxc, maxc) = T(1);
1021 result.set_row(maxc, result.row(maxc) * rpivot);
1022
1023 // reduce rows (except pivot):
1024 for (size_t j = 0; j < N; ++j) {
1025 if (j != maxc) {
1026 T dum = result(j, maxc);
1027 result(j, maxc) = T(0);
1028 for (size_t k = 0; k < N; ++k)
1029 result(j, k) -= result(maxc, k) * dum;
1030 }
1031 }
1032 }
1033
1034 // undo column swap:
1035 size_t i = N;
1036 do {
1037 --i;
1038 if (indxr[i] != indxc[i]) // if swap occurred
1039 result.swap_cols(indxr[i], indxc[i]);
1040 } while (i != 0);
1041
1042 return result;
1043 }
1044
1046 template <typename T>
1047 inline Mat2<T> inverse(const Mat2<T> &m) {
1048 Mat2<T> result;
1049 result(0, 0) = m(1, 1);
1050 result(0, 1) = -m(0, 1);
1051 result(1, 0) = -m(1, 0);
1052 result(1, 1) = m(0, 0);
1053
1054 T det = determinant(m);
1055 det = T(1) / det;
1056 result *= det;
1057
1058 return result;
1059 }
1060
1062 template <typename T>
1063 inline Mat3<T> inverse(const Mat3<T> &m) {
1064 Mat3<T> result;
1065 result(0, 0) = (m(1, 1) * m(2, 2) - m(2, 1) * m(1, 2));
1066 result(0, 1) = -(m(0, 1) * m(2, 2) - m(0, 2) * m(2, 1));
1067 result(0, 2) = (m(0, 1) * m(1, 2) - m(0, 2) * m(1, 1));
1068 result(1, 0) = -(m(1, 0) * m(2, 2) - m(1, 2) * m(2, 0));
1069 result(1, 1) = (m(0, 0) * m(2, 2) - m(0, 2) * m(2, 0));
1070 result(1, 2) = -(m(0, 0) * m(1, 2) - m(1, 0) * m(0, 2));
1071 result(2, 0) = (m(1, 0) * m(2, 1) - m(2, 0) * m(1, 1));
1072 result(2, 1) = -(m(0, 0) * m(2, 1) - m(2, 0) * m(0, 1));
1073 result(2, 2) = (m(0, 0) * m(1, 1) - m(1, 0) * m(0, 1));
1074
1075 T det = determinant(m);
1076 det = T(1) / det;
1077 result *= det;
1078
1079 return result;
1080 }
1081
1083 template <typename T>
1084 inline Mat4<T> inverse(const Mat4<T> &m) {
1085 Mat4<T> result;
1086 result(0, 0) = m(1, 2) * m(2, 3) * m(3, 1) - m(1, 3) * m(2, 2) * m(3, 1) + m(1, 3) * m(2, 1) * m(3, 2) - m(1, 1) * m(2, 3) * m(3, 2) - m(1, 2) * m(2, 1) * m(3, 3) + m(1, 1) * m(2, 2) * m(3, 3);
1087 result(0, 1) = m(0, 3) * m(2, 2) * m(3, 1) - m(0, 2) * m(2, 3) * m(3, 1) - m(0, 3) * m(2, 1) * m(3, 2) + m(0, 1) * m(2, 3) * m(3, 2) + m(0, 2) * m(2, 1) * m(3, 3) - m(0, 1) * m(2, 2) * m(3, 3);
1088 result(0, 2) = m(0, 2) * m(1, 3) * m(3, 1) - m(0, 3) * m(1, 2) * m(3, 1) + m(0, 3) * m(1, 1) * m(3, 2) - m(0, 1) * m(1, 3) * m(3, 2) - m(0, 2) * m(1, 1) * m(3, 3) + m(0, 1) * m(1, 2) * m(3, 3);
1089 result(0, 3) = m(0, 3) * m(1, 2) * m(2, 1) - m(0, 2) * m(1, 3) * m(2, 1) - m(0, 3) * m(1, 1) * m(2, 2) + m(0, 1) * m(1, 3) * m(2, 2) + m(0, 2) * m(1, 1) * m(2, 3) - m(0, 1) * m(1, 2) * m(2, 3);
1090 result(1, 0) = m(1, 3) * m(2, 2) * m(3, 0) - m(1, 2) * m(2, 3) * m(3, 0) - m(1, 3) * m(2, 0) * m(3, 2) + m(1, 0) * m(2, 3) * m(3, 2) + m(1, 2) * m(2, 0) * m(3, 3) - m(1, 0) * m(2, 2) * m(3, 3);
1091 result(1, 1) = m(0, 2) * m(2, 3) * m(3, 0) - m(0, 3) * m(2, 2) * m(3, 0) + m(0, 3) * m(2, 0) * m(3, 2) - m(0, 0) * m(2, 3) * m(3, 2) - m(0, 2) * m(2, 0) * m(3, 3) + m(0, 0) * m(2, 2) * m(3, 3);
1092 result(1, 2) = m(0, 3) * m(1, 2) * m(3, 0) - m(0, 2) * m(1, 3) * m(3, 0) - m(0, 3) * m(1, 0) * m(3, 2) + m(0, 0) * m(1, 3) * m(3, 2) + m(0, 2) * m(1, 0) * m(3, 3) - m(0, 0) * m(1, 2) * m(3, 3);
1093 result(1, 3) = m(0, 2) * m(1, 3) * m(2, 0) - m(0, 3) * m(1, 2) * m(2, 0) + m(0, 3) * m(1, 0) * m(2, 2) - m(0, 0) * m(1, 3) * m(2, 2) - m(0, 2) * m(1, 0) * m(2, 3) + m(0, 0) * m(1, 2) * m(2, 3);
1094 result(2, 0) = m(1, 1) * m(2, 3) * m(3, 0) - m(1, 3) * m(2, 1) * m(3, 0) + m(1, 3) * m(2, 0) * m(3, 1) - m(1, 0) * m(2, 3) * m(3, 1) - m(1, 1) * m(2, 0) * m(3, 3) + m(1, 0) * m(2, 1) * m(3, 3);
1095 result(2, 1) = m(0, 3) * m(2, 1) * m(3, 0) - m(0, 1) * m(2, 3) * m(3, 0) - m(0, 3) * m(2, 0) * m(3, 1) + m(0, 0) * m(2, 3) * m(3, 1) + m(0, 1) * m(2, 0) * m(3, 3) - m(0, 0) * m(2, 1) * m(3, 3);
1096 result(2, 2) = m(0, 1) * m(1, 3) * m(3, 0) - m(0, 3) * m(1, 1) * m(3, 0) + m(0, 3) * m(1, 0) * m(3, 1) - m(0, 0) * m(1, 3) * m(3, 1) - m(0, 1) * m(1, 0) * m(3, 3) + m(0, 0) * m(1, 1) * m(3, 3);
1097 result(2, 3) = m(0, 3) * m(1, 1) * m(2, 0) - m(0, 1) * m(1, 3) * m(2, 0) - m(0, 3) * m(1, 0) * m(2, 1) + m(0, 0) * m(1, 3) * m(2, 1) + m(0, 1) * m(1, 0) * m(2, 3) - m(0, 0) * m(1, 1) * m(2, 3);
1098 result(3, 0) = m(1, 2) * m(2, 1) * m(3, 0) - m(1, 1) * m(2, 2) * m(3, 0) - m(1, 2) * m(2, 0) * m(3, 1) + m(1, 0) * m(2, 2) * m(3, 1) + m(1, 1) * m(2, 0) * m(3, 2) - m(1, 0) * m(2, 1) * m(3, 2);
1099 result(3, 1) = m(0, 1) * m(2, 2) * m(3, 0) - m(0, 2) * m(2, 1) * m(3, 0) + m(0, 2) * m(2, 0) * m(3, 1) - m(0, 0) * m(2, 2) * m(3, 1) - m(0, 1) * m(2, 0) * m(3, 2) + m(0, 0) * m(2, 1) * m(3, 2);
1100 result(3, 2) = m(0, 2) * m(1, 1) * m(3, 0) - m(0, 1) * m(1, 2) * m(3, 0) - m(0, 2) * m(1, 0) * m(3, 1) + m(0, 0) * m(1, 2) * m(3, 1) + m(0, 1) * m(1, 0) * m(3, 2) - m(0, 0) * m(1, 1) * m(3, 2);
1101 result(3, 3) = m(0, 1) * m(1, 2) * m(2, 0) - m(0, 2) * m(1, 1) * m(2, 0) + m(0, 2) * m(1, 0) * m(2, 1) - m(0, 0) * m(1, 2) * m(2, 1) - m(0, 1) * m(1, 0) * m(2, 2) + m(0, 0) * m(1, 1) * m(2, 2);
1102
1103 T det = determinant(m);
1104 det = T(1) / det;
1105 result *= det;
1106
1107 return result;
1108 }
1109
1110 /*----------------------------------------------------------------------------*/
1111 template <size_t M, size_t N, typename T>
1112 inline Mat<N, M, T> tensor(const Vec<M, T> &u, const Vec<N, T> &v) {
1113 Mat<N, 1, T> mu(u); // column vector u
1114 Mat<1, M, T> mv(v); // row vector v
1115 return mu * mv; // use matrix multiplication
1116 }
1117
1118 /*----------------------------------------------------------------------------*/
1119 /* Adapted from "Numerical Recipes in C" (Press et al.) */
1120 template <size_t N, size_t M, typename T>
1122 const Mat<N, N, T> &a,
1123 const Mat<N, M, T> &b,
1124 Mat<N, N, T> *ainv,
1125 Mat<N, M, T> *x) {
1126
1127 size_t indxc[N], indxr[N], ipiv[N];
1128
1129 Mat<N, N, T> &amat = *ainv;
1130 amat = a; // copy A
1131 Mat<N, M, T> &bmat = *x;
1132 bmat = b; // copy B
1133
1134 for (size_t i = 0; i < N; ++i)
1135 ipiv[i] = 0;
1136
1137 for (size_t i = 0; i < N; ++i) { // for each column
1138 T max = T(0);
1139 size_t maxc, maxr;
1140 for (size_t j = 0; j < N; ++j) { // search for pivot
1141 if (ipiv[j] != 1) {
1142 for (size_t k = 0; k < N; ++k) { // for each column
1143 if (ipiv[k] == 0) {
1144 T element = std::abs(amat(j, k));
1145 if (element > max) {
1146 max = element;
1147 maxr = j;
1148 maxc = k;
1149 }
1150 }
1151 }
1152 }
1153 }
1154 ++ipiv[maxc];
1155
1156 if (maxr != maxc) {
1157 amat.swap_rows(maxr, maxc);
1158 bmat.swap_rows(maxr, maxc);
1159 }
1160
1161 indxr[i] = maxr;
1162 indxc[i] = maxc;
1163
1164 // check for singular matrix:
1165 if (std::abs(amat(maxc, maxc)) < epsilon<T>())
1166 return false;
1167
1168 // multiply row by 1/pivot:
1169 T rpivot = T(1) / amat(maxc, maxc);
1170 amat(maxc, maxc) = T(1);
1171 amat.set_row(maxc, amat.row(maxc) * rpivot);
1172 bmat.set_row(maxc, bmat.row(maxc) * rpivot);
1173
1174 // reduce rows (except pivot):
1175 for (size_t j = 0; j < N; ++j) {
1176 if (j != maxc) {
1177 T dum = amat(j, maxc);
1178 amat(j, maxc) = T(0);
1179 for (size_t k = 0; k < N; ++k)
1180 amat(j, k) -= amat(maxc, k) * dum;
1181 for (size_t k = 0; k < M; ++k)
1182 bmat(j, k) -= bmat(maxc, k) * dum;
1183 }
1184 }
1185 }
1186
1187 // unscramble a by interchanging columns:
1188 size_t i = N;
1189 do {
1190 --i;
1191 if (indxr[i] != indxc[i]) // if swap occurred
1192 amat.swap_cols(indxr[i], indxc[i]);
1193 } while (i != 0);
1194
1195 return true;
1196 }
1197
1198 /*----------------------------------------------------------------------------*/
1199 /* Adapted from "Numerical Recipes in C" (Press et al.) */
1200 template <size_t N, typename T>
1201 inline bool lu_decomposition(
1202 const Mat<N, N, T> &a,
1203 Mat<N, N, T> *alu,
1204 Vec<N, T> *rowp,
1205 T *d) {
1206
1207 Vec<N, T> scalev; // stores implicit scaling of each row
1208 *d = T(1); // no rows changed
1209 Mat<N, N, T> &amat = *alu;
1210 amat = a; // copy a
1211
1212 for (size_t i = 0; i < N; ++i) { // for each row
1213 // get implicit scaling:
1214 T max = T(0);
1215 for (size_t j = 0; j < N; ++j) {
1216 T element = std::abs(amat(i, j));
1217 if (element > max)
1218 max = element;
1219 }
1220
1221 // check for singular matrix:
1222 if (std::abs(max) < min<T>())
1223 return false;
1224
1225 scalev[i] = T(1) / max; // save scaling factor
1226 }
1227
1228 for (size_t j = 0; j < N; ++j) { // for each column
1229 for (size_t i = 0; i < j; ++i) {
1230 T sum = amat(i, j);
1231 for (size_t k = 0; k < i; ++k)
1232 sum -= amat(i, k) * amat(k, j);
1233 amat(i, j) = sum;
1234 }
1235
1236 T max = T(0);
1237 size_t imax;
1238 for (size_t i = j; i < N; ++i) {
1239 T sum = amat(i, j);
1240 for (size_t k = 0; k < j; ++k)
1241 sum -= amat(i, k) * amat(k, j);
1242 amat(i, j) = sum;
1243
1244 T dum = scalev[i] * std::abs(sum);
1245 if (dum >= max) {
1246 max = dum;
1247 imax = i;
1248 }
1249 }
1250
1251 if (j != imax) { // if we need to swap rows
1252 amat.swap_rows(imax, j);
1253 scalev[imax] = scalev[j]; // also swap scale factor
1254 *d = -(*d); // change parity of d
1255 }
1256 (*rowp)[j] = (T)imax;
1257
1258 // check for singular matrix:
1259 if (std::abs(amat(j, j)) < epsilon<T>())
1260 return false;
1261
1262 // divide by the pivot element:
1263 if (j != N) {
1264 T dum = T(1) / amat(j, j);
1265 for (size_t i = j + 1; i < N; ++i)
1266 amat(i, j) *= dum;
1267 }
1268
1269 }
1270
1271 return true;
1272 }
1273
1274 /*----------------------------------------------------------------------------*/
1275 /* Adapted from "Numerical Recipes in C" (Press et al.) */
1276 template <size_t N, typename T>
1278 const Mat<N, N, T> &alu,
1279 const Vec<N, T> &rowp,
1280 const Vec<N, T> &b,
1281 Vec<N, T> *x
1282 )
1283 {
1284 Vec<N, T> &result = *x;
1285 result = b; // copy b to result
1286
1287 size_t ii = 0;
1288 for (size_t i = 0; i < N; ++i) {
1289 size_t ip = (size_t)rowp[i];
1290 assert(ip < N);
1291
1292 T sum = result[ip];
1293 result[ip] = result[i];
1294 if (ii != 0) {
1295 for (size_t j = ii - 1; j < i; ++j)
1296 sum -= alu(i, j) * result[j];
1297 }
1298 else if (std::abs(sum) > epsilon<T>()) {
1299 ii = i + 1;
1300 }
1301 result[i] = sum;
1302 }
1303
1304 size_t i = N;
1305 do {
1306 --i;
1307 T sum = result[i];
1308 for (size_t j = i + 1; j < N; ++j)
1309 sum -= alu(i, j) * result[j];
1310 result[i] = sum / alu(i, i);
1311 } while (i != 0);
1312 }
1313
1314
1315 /*----------------------------------------------------------------------------*/
1316 // Adapted from Template Numerical Toolkit.
1317 template<size_t N, typename FT>
1319 bool spd = true;
1320 for (size_t j = 0; j < N; ++j) {
1321 FT d = 0;
1322 for (size_t k = 0; k < j; ++k) {
1323 FT s = 0;
1324 for (size_t i = 0; i < k; ++i)
1325 s += L(k, i) * L(j, i);
1326
1327 L(j, k) = s = (A(j, k) - s) / L(k, k);
1328 d = d + s * s;
1329 spd = spd && (A(k, j) == A(j, k));
1330 }
1331
1332 d = A(j, j) - d;
1333 spd = spd && (d > 0);
1334
1335 L(j, j) = std::sqrt(d > 0 ? d : 0);
1336 for (size_t k = j + 1; k < N; ++k)
1337 L(j, k) = 0;
1338 }
1339 return spd;
1340 }
1341
1342 /*----------------------------------------------------------------------------*/
1343 // Adapted from Template Numerical Toolkit.
1344 template <size_t N, typename FT>
1345 void cholesky_solve(const Mat<N, N, FT> &L, const Vec<N, FT> &b, Vec<N, FT>& x) {
1346 x = b;
1347 // solve L*y = b
1348 for (size_t k = 0; k < N; ++k) {
1349 for (size_t i = 0; i < k; ++i)
1350 x[k] -= x[i] * L(k, i);
1351 x[k] /= L(k, k);
1352 }
1353
1354 // solve L'*x = y
1355 for (int k = N - 1; k >= 0; --k) { // attention: size_t will not work
1356 for (size_t i = k + 1; i < N; ++i)
1357 x[k] -= x[i] * L(i, k);
1358 x[k] /= L(k, k);
1359 }
1360 }
1361
1362 template<size_t N, size_t M, typename T>
1364 X = B;
1365 // solve L_*Y = B
1366 for (size_t j = 0; j < M; ++j) {
1367 for (size_t k = 0; k < N; ++k) {
1368 for (size_t i = 0; i < k; ++i)
1369 X(k, j) -= X(i, j) * L(k, i);
1370
1371 X(k, j) /= L(k, k);
1372 }
1373 }
1374
1375 // solve L'*X = Y
1376 for (size_t j = 0; j < M; ++j) {
1377 for (int k = N - 1; k >= 0; --k) { // attention: size_t will not work
1378 for (size_t i = k + 1; i < N; ++i)
1379 X(k, j) -= X(i, j) * L(i, k);
1380
1381 X(k, j) /= L(k, k);
1382 }
1383 }
1384 }
1385
1386 /*----------------------------------------------------------------------------*/
1387
1389 template <size_t N, size_t M, typename T> inline
1390 std::ostream& operator<< (std::ostream& output, const Mat<N, M, T>& m) {
1391 output << std::fixed << std::setprecision(8);
1392 const char sep = ' ';
1393 for (int i = 0; i < N; i++) {
1394 for (int j = 0; j < M; j++) {
1395 output << sep << std::setw(7) << m(i, j);
1396 }
1397 output << std::endl;
1398 }
1399 output << std::resetiosflags(std::ios_base::fixed | std::ios_base::floatfield);
1400 return output;
1401 }
1402
1403 /*----------------------------------------------------------------------------*/
1405 template <size_t N, size_t M, typename T>
1406 std::istream& operator>> (std::istream& input, Mat<N, M, T>& m) {
1407 for (int i = 0; i < N; i++) {
1408 for (int j = 0; j < M; j++) {
1409 input >> m(i, j);
1410 }
1411 }
1412 return input;
1413 }
1414
1415
1416 /*----------------------------------------------------------------------------*/
1417
1419 template<size_t N, typename FT>
1421 return Mat<N, 1, FT>(v.data());
1422 }
1423
1424
1426 template<size_t N, typename FT>
1428 return Mat<1, N, FT>(v.data());
1429 }
1430
1431 /*----------------------------------------------------------------------------*/
1432
1434 template <size_t N, size_t M, typename T>
1435 inline bool has_nan(const Mat<N, M, T>& m) {
1436 for (int i = 0; i < N; i++) {
1437 for (int j = 0; j < M; j++) {
1438 if (std::isnan(m(i, j)) || std::isinf(m(i, j)))
1439 return true;
1440 }
1441 }
1442 return false;
1443 }
1444
1445
1446 /*******************************************************************************
1447
1448 -------------------------------- 2x2 matrix ----------------------------------
1449
1450 *******************************************************************************/
1451
1458 template <typename T>
1459 class Mat2 : public Mat<2, 2, T>
1460 {
1461 public:
1468 Mat2() = default;
1469
1471 explicit Mat2(T s);
1472
1474 Mat2(const Mat<2, 2, T> &rhs);
1475
1480 explicit Mat2(const Mat<3, 3, T> &rhs);
1481
1486 Mat2(
1487 T s00, T s01,
1488 T s10, T s11
1489 );
1490
1492 explicit Mat2(const T *m);
1493
1498 Mat2(
1499 const Vec<2, T> &x,
1500 const Vec<2, T> &y
1501 );
1502
1508 static Mat2<T> rotation(T angle);
1509
1514 static Mat2<T> scale(T s);
1515
1521 static Mat2<T> scale(T x, T y);
1522
1523 }; // class Mat2
1524
1525 /*----------------------------------------------------------------------------*/
1526 template <typename T>
1527 inline Mat2<T>::Mat2(T s) {
1528 for (size_t i = 0; i < 2; ++i) {
1529 for (size_t j = 0; j < 2; ++j)
1530 if (i == j)
1531 (*this)(i, j) = T(s);
1532 else
1533 (*this)(i, j) = T(0);
1534 }
1535 }
1536
1537 /*----------------------------------------------------------------------------*/
1538 template <typename T>
1539 inline Mat2<T>::Mat2(const Mat<2, 2, T> &rhs) {
1540 for (size_t i = 0; i < 4; ++i)
1541 (*this).m_[i] = rhs[i];
1542 }
1543
1544 /*----------------------------------------------------------------------------*/
1545 template <typename T>
1546 inline Mat2<T>::Mat2(const Mat<3, 3, T> &rhs) {
1547 (*this)(0, 0) = rhs(0, 0); (*this)(0, 1) = rhs(0, 1);
1548 (*this)(1, 0) = rhs(1, 0); (*this)(1, 1) = rhs(1, 1);
1549 }
1550
1551 /*----------------------------------------------------------------------------*/
1552 template <typename T>
1554 T s00, T s01,
1555 T s10, T s11
1556 ) {
1557 (*this)(0, 0) = s00; (*this)(0, 1) = s01;
1558 (*this)(1, 0) = s10; (*this)(1, 1) = s11;
1559 }
1560
1561 /*----------------------------------------------------------------------------*/
1562 template <typename T>
1563 inline Mat2<T>::Mat2(const T *m) {
1564 assert(m != 0);
1565
1566 #ifdef MATRIX_ROW_MAJOR
1567 (*this)(0, 0) = m[0]; (*this)(0, 1) = m[1];
1568 (*this)(1, 0) = m[2]; (*this)(1, 1) = m[3];
1569 #else
1570 (*this)(0, 0) = m[0]; (*this)(0, 1) = m[2];
1571 (*this)(1, 0) = m[1]; (*this)(1, 1) = m[3];
1572 #endif
1573 }
1574
1575 /*----------------------------------------------------------------------------*/
1576 template <typename T>
1577 inline Mat2<T>::Mat2(const Vec<2, T> &x, const Vec<2, T> &y) {
1578 #ifdef MATRIX_ROW_MAJOR
1579 (*this).set_row(0, x);
1580 (*this).set_row(1, y);
1581 #else
1582 (*this).set_col(0, x);
1583 (*this).set_col(1, y);
1584 #endif
1585 }
1586
1587 /*----------------------------------------------------------------------------*/
1588 template <typename T>
1589 inline Mat2<T> Mat2<T>::rotation(T angle) {
1590 return Mat2<T>(
1591 std::cos(angle), -std::sin(angle),
1592 std::sin(angle), std::cos(angle)
1593 );
1594 }
1595
1596 /*----------------------------------------------------------------------------*/
1597 template <typename T>
1599 return Mat2<T>(
1600 s, T(0),
1601 T(0), s
1602 );
1603 }
1604
1605 /*----------------------------------------------------------------------------*/
1606 template <typename T>
1607 inline Mat2<T> Mat2<T>::scale(T x, T y) {
1608 return Mat2<T>(
1609 x, T(0),
1610 T(0), y
1611 );
1612 }
1613
1614 /*******************************************************************************
1615
1616 -------------------------------- 3x3 matrix ----------------------------------
1617
1618 *******************************************************************************/
1619
1620 template <typename T> class Quat;
1621
1627 template <typename T>
1628 class Mat3 : public Mat<3, 3, T> {
1629 public:
1636 Mat3() = default;
1637
1639 explicit Mat3(T s);
1640
1642 Mat3(const Mat<3, 3, T> &rhs);
1643
1648 explicit Mat3(const Mat<4, 4, T> &rhs);
1649
1654 Mat3(
1655 T s00, T s01, T s02,
1656 T s10, T s11, T s12,
1657 T s20, T s21, T s22
1658 );
1659
1661 explicit Mat3(const T *m);
1662
1667 Mat3(
1668 const Vec<3, T> &x,
1669 const Vec<3, T> &y,
1670 const Vec<3, T> &z
1671 );
1672
1677 explicit Mat3(const Mat<2, 2, T> &rhs);
1678
1680 Mat2<T> sub() const;
1681
1686 static Mat3<T> scale(T s);
1687
1694 static Mat3<T> scale(T x, T y, T z);
1695
1704 static Mat3<T> rotation(const Vec<3, T> &axis, T angle);
1705
1714 static Mat3<T> rotation(const Vec<3, T> &axis_angle);
1715
1719 static Mat3<T> rotation(const Quat<T> &q);
1720
1728 static Mat3<T> rotation(T x, T y, T z, int order = 312);
1729
1730 }; // class Mat3
1731
1732
1733 /*******************************************************************************
1734
1735 IMPLEMENTATION:
1736
1737 *******************************************************************************/
1738
1739 /*----------------------------------------------------------------------------*/
1740 template <typename T>
1741 inline Mat3<T>::Mat3(T s) {
1742 for (size_t i = 0; i < 3; ++i) {
1743 for (size_t j = 0; j < 3; ++j)
1744 if (i == j)
1745 (*this)(i, j) = T(s);
1746 else
1747 (*this)(i, j) = T(0);
1748 }
1749 }
1750
1751 /*----------------------------------------------------------------------------*/
1752 template <typename T>
1753 inline Mat3<T>::Mat3(const Mat<3, 3, T> &rhs) {
1754 for (size_t i = 0; i < 9; ++i)
1755 (*this).m_[i] = rhs[i];
1756 }
1757
1758 /*----------------------------------------------------------------------------*/
1759 template <typename T>
1760 inline Mat3<T>::Mat3(const Mat<4, 4, T> &rhs) {
1761 (*this)(0, 0) = rhs(0, 0); (*this)(0, 1) = rhs(0, 1); (*this)(0, 2) = rhs(0, 2);
1762 (*this)(1, 0) = rhs(1, 0); (*this)(1, 1) = rhs(1, 1); (*this)(1, 2) = rhs(1, 2);
1763 (*this)(2, 0) = rhs(2, 0); (*this)(2, 1) = rhs(2, 1); (*this)(2, 2) = rhs(2, 2);
1764 }
1765
1766 /*----------------------------------------------------------------------------*/
1767 template <typename T>
1769 T s00, T s01, T s02,
1770 T s10, T s11, T s12,
1771 T s20, T s21, T s22
1772 ) {
1773 (*this)(0, 0) = s00; (*this)(0, 1) = s01; (*this)(0, 2) = s02;
1774 (*this)(1, 0) = s10; (*this)(1, 1) = s11; (*this)(1, 2) = s12;
1775 (*this)(2, 0) = s20; (*this)(2, 1) = s21; (*this)(2, 2) = s22;
1776 }
1777
1778 /*----------------------------------------------------------------------------*/
1779 template <typename T>
1780 inline Mat3<T>::Mat3(const T *m) {
1781 assert(m != 0);
1782
1783 #ifdef MATRIX_ROW_MAJOR
1784 (*this)(0, 0) = m[0]; (*this)(0, 1) = m[1]; (*this)(0, 2) = m[2];
1785 (*this)(1, 0) = m[3]; (*this)(1, 1) = m[4]; (*this)(1, 2) = m[5];
1786 (*this)(2, 0) = m[6]; (*this)(2, 1) = m[7]; (*this)(2, 2) = m[8];
1787 #else
1788 (*this)(0, 0) = m[0]; (*this)(0, 1) = m[3]; (*this)(0, 2) = m[6];
1789 (*this)(1, 0) = m[1]; (*this)(1, 1) = m[4]; (*this)(1, 2) = m[7];
1790 (*this)(2, 0) = m[2]; (*this)(2, 1) = m[5]; (*this)(2, 2) = m[8];
1791 #endif
1792 }
1793
1794 /*----------------------------------------------------------------------------*/
1795 template <typename T>
1796 inline Mat3<T>::Mat3(const Vec<3, T> &x, const Vec<3, T> &y, const Vec<3, T> &z) {
1797 #ifdef MATRIX_ROW_MAJOR
1798 (*this).set_row(0, x);
1799 (*this).set_row(1, y);
1800 (*this).set_row(2, z);
1801 #else
1802 (*this).set_col(0, x);
1803 (*this).set_col(1, y);
1804 (*this).set_col(2, z);
1805 #endif
1806 }
1807
1808 /*----------------------------------------------------------------------------*/
1809 template <typename T>
1810 inline Mat3<T>::Mat3(const Mat<2, 2, T> &rhs) {
1811 (*this)(0, 0) = rhs(0, 1); (*this)(0, 0) = rhs(0, 1); (*this)(0, 2) = T(0);
1812 (*this)(1, 0) = rhs(1, 1); (*this)(1, 0) = rhs(1, 1); (*this)(1, 2) = T(0);
1813 (*this)(2, 0) = T(0); (*this)(2, 0) = T(0); (*this)(2, 2) = T(1);
1814 }
1815
1816 /*----------------------------------------------------------------------------*/
1817 template <typename T>
1818 inline Mat2<T> Mat3<T>::sub() const {
1819 Mat2<T> mat;
1820 for (size_t i = 0; i < 2; i++) {
1821 for (size_t j = 0; j < 2; j++) {
1822 mat(i, j) = (*this)(i, j);
1823 }
1824 }
1825 return mat;
1826 }
1827
1828 /*----------------------------------------------------------------------------*/
1829 template <typename T>
1831 return Mat3<T>(
1832 s, T(0), T(0),
1833 T(0), s, T(0),
1834 T(0), T(0), s
1835 );
1836 }
1837
1838 /*----------------------------------------------------------------------------*/
1839 template <typename T>
1840 inline Mat3<T> Mat3<T>::scale(T x, T y, T z) {
1841 return Mat3<T>(
1842 x, T(0), T(0),
1843 T(0), y, T(0),
1844 T(0), T(0), z
1845 );
1846 }
1847
1848 /*----------------------------------------------------------------------------*/
1849 template <typename T>
1850 inline Mat3<T> Mat3<T>::rotation(const Vec<3, T> &axis, T angle) {
1851 assert(std::abs(axis.length() - 1) < epsilon<T>());
1852
1853 // cross-product matrix of axis:
1854 const Mat3<T> cpm(
1855 T(0), -axis[2], axis[1],
1856 axis[2], T(0), -axis[0],
1857 -axis[1], axis[0], T(0)
1858 );
1859
1860 // axis-axis tensor product:
1861 const Mat3<T> tpm = tensor(axis, axis);
1862
1863 // trig coefficients:
1864 const T c = std::cos(angle);
1865 const T rc = T(1) - c;
1866 const T s = std::sin(angle);
1867
1868 return Mat3<T>::identity() * c + cpm * s + tpm * rc;
1869 }
1870
1871 /*----------------------------------------------------------------------------*/
1872 template <typename T>
1873 inline Mat3<T> Mat3<T>::rotation(const Vec<3, T> &axis_angle) {
1874 const T len = axis_angle.length();
1875 return rotation(axis_angle/len, len);
1876 }
1877
1878 /*----------------------------------------------------------------------------*/
1879 template <typename T>
1881 // input must be unit quaternion
1882 assert(std::abs(q.length() - 1) < epsilon<T>());
1883 const T x = q.x;
1884 const T y = q.y;
1885 const T z = q.z;
1886 const T w = q.w;
1887 Mat3<T> m;
1888 m(0, 0) = 1.0 - 2.0 * (y*y + z*z); m(0, 1) = 2.0 * (x*y - w*z); m(0, 2) = 2.0 * (x*z + w*y);
1889 m(1, 0) = 2.0 * (x*y + w*z); m(1, 1) = 1.0 - 2.0 * (x*x + z*z); m(1, 2) = 2.0 * (y*z - w*x);
1890 m(2, 0) = 2.0 * (x*z - w*y); m(2, 1) = 2.0 * (y*z + w*x); m(2, 2) = 1.0 - 2.0 * (x*x + y*y);
1891 return m;
1892 }
1893
1894 /*----------------------------------------------------------------------------*/
1895 template <typename T>
1896 inline Mat3<T> Mat3<T>::rotation(T x, T y, T z, int order) {
1897 // The code is not optimized. The final matrix can be directly derived.
1898 // http://www.songho.ca/opengl/gl_anglestoaxes.html
1899 // http://inside.mines.edu/fs_home/gmurray/ArbitraryAxisRotation/
1900 Mat3<T> rx; // Rotation about X-axis (Pitch)
1901 rx(0, 0) = T(1); rx(0, 1) = T(0); rx(0, 2) = T(0);
1902 rx(1, 0) = T(0); rx(1, 1) = std::cos(x); rx(1, 2) = -std::sin(x);
1903 rx(2, 0) = T(0); rx(2, 1) = std::sin(x); rx(2, 2) = std::cos(x);
1904
1905 Mat3<T> ry; // Rotation about Y-axis (Yaw, Heading)
1906 ry(0, 0) = std::cos(y); ry(0, 1) = T(0); ry(0, 2) = std::sin(y);
1907 ry(1, 0) = T(0); ry(1, 1) = T(1); ry(1, 2) = T(0);
1908 ry(2, 0) = -std::sin(y); ry(2, 1) = T(0); ry(2, 2) = std::cos(y);
1909
1910 Mat3<T> rz; // Rotation about Z-axis (Roll)
1911 rz(0, 0) = std::cos(z); rz(0, 1) = -std::sin(z); rz(0, 2) = T(0);
1912 rz(1, 0) = std::sin(z); rz(1, 1) = std::cos(z); rz(1, 2) = T(0);
1913 rz(2, 0) = T(0); rz(2, 1) = T(0); rz(2, 2) = T(1);
1914
1915 switch (order) {
1916 case 123: return rz * ry * rx;
1917 case 132: return ry * rz * rx;
1918 case 213: return rz * rx * ry;
1919 case 231: return rx * rz * ry;
1920 case 312: return ry * rx * rz;
1921 case 321: return rx * ry * rz;
1922 default:
1923 LOG(ERROR) << "invalid rotation order";
1924 return rx * rz * ry;
1925 }
1926 }
1927
1928 /*******************************************************************************
1929
1930 -------------------------------- 4x4 matrix ----------------------------------
1931
1932 *******************************************************************************/
1933
1939 template <typename T>
1940 class Mat4 : public Mat<4, 4, T> {
1941 public:
1948 Mat4() = default;
1949
1951 explicit Mat4(T s);
1952
1954 Mat4(const Mat<4, 4, T> &rhs);
1955
1961 T s00, T s01, T s02, T s03,
1962 T s10, T s11, T s12, T s13,
1963 T s20, T s21, T s22, T s23,
1964 T s30, T s31, T s32, T s33
1965 );
1966
1968 explicit Mat4(const T *m);
1969
1975 const Vec<4, T> &x,
1976 const Vec<4, T> &y,
1977 const Vec<4, T> &z,
1978 const Vec<4, T> &w
1979 );
1980
1985 explicit Mat4(const Mat<3, 3, T> &rhs);
1986
1988 Mat4(const Vec<3, T> &s, const Quat<T> &r, const Vec<3, T> &t);
1989
1990 // the upper-left 3x3 sub-matrix,
1991 Mat3<T> sub() const;
1992
1997 static Mat4<T> scale(T s);
1998
2006 static Mat4<T> scale(const Vec<4, T>& s); // set w to 1 for 3D scaling
2007 static Mat4<T> scale(T x, T y, T z, T w); // set w to 1 for 3D scaling
2008
2017 static Mat4<T> rotation(const Vec<3, T> &axis, T angle);
2018
2027 static Mat4<T> rotation(const Vec<3, T> &axis_angle);
2028
2033 static Mat4<T> rotation(const Quat<T> &q);
2034
2043 static Mat4<T> rotation(T x, T y, T z, int order = 312);
2044
2050 static Mat4<T> translation(T x, T y, T z);
2051
2052 }; // class Mat4
2053
2054
2055 /*******************************************************************************
2056
2057 IMPLEMENTATION:
2058
2059 *******************************************************************************/
2060
2061 /*----------------------------------------------------------------------------*/
2062 template <typename T>
2063 inline Mat4<T>::Mat4(T s) {
2064 for (size_t i = 0; i < 4; ++i) {
2065 for (size_t j = 0; j < 4; ++j)
2066 if (i == j)
2067 (*this)(i, j) = T(s);
2068 else
2069 (*this)(i, j) = T(0);
2070 }
2071 }
2072
2073 /*----------------------------------------------------------------------------*/
2074 template <typename T>
2075 inline Mat4<T>::Mat4(const Mat<4, 4, T> &rhs) {
2076 for (size_t i = 0; i < 16; ++i)
2077 (*this).m_[i] = rhs[i];
2078 }
2079
2080 /*----------------------------------------------------------------------------*/
2081 template <typename T>
2083 T s00, T s01, T s02, T s03,
2084 T s10, T s11, T s12, T s13,
2085 T s20, T s21, T s22, T s23,
2086 T s30, T s31, T s32, T s33
2087 ) {
2088 (*this)(0, 0) = s00; (*this)(0, 1) = s01; (*this)(0, 2) = s02; (*this)(0, 3) = s03;
2089 (*this)(1, 0) = s10; (*this)(1, 1) = s11; (*this)(1, 2) = s12; (*this)(1, 3) = s13;
2090 (*this)(2, 0) = s20; (*this)(2, 1) = s21; (*this)(2, 2) = s22; (*this)(2, 3) = s23;
2091 (*this)(3, 0) = s30; (*this)(3, 1) = s31; (*this)(3, 2) = s32; (*this)(3, 3) = s33;
2092 }
2093
2094 /*----------------------------------------------------------------------------*/
2095 template <typename T>
2096 inline Mat4<T>::Mat4(const T *m) {
2097 assert(m != 0);
2098
2099 #ifdef MATRIX_ROW_MAJOR
2100 (*this)(0, 0) = m[0]; (*this)(0, 1) = m[1]; (*this)(0, 2) = m[2]; (*this)(0, 3) = m[3];
2101 (*this)(1, 0) = m[4]; (*this)(1, 1) = m[5]; (*this)(1, 2) = m[6]; (*this)(1, 3) = m[7];
2102 (*this)(2, 0) = m[8]; (*this)(2, 1) = m[9]; (*this)(2, 2) = m[10]; (*this)(2, 3) = m[11];
2103 (*this)(3, 0) = m[12]; (*this)(3, 1) = m[13]; (*this)(3, 2) = m[14]; (*this)(3, 3) = m[15];
2104 #else
2105 (*this)(0, 0) = m[0]; (*this)(0, 1) = m[4]; (*this)(0, 2) = m[8]; (*this)(0, 3) = m[12];
2106 (*this)(1, 0) = m[1]; (*this)(1, 1) = m[5]; (*this)(1, 2) = m[9]; (*this)(1, 3) = m[13];
2107 (*this)(2, 0) = m[2]; (*this)(2, 1) = m[6]; (*this)(2, 2) = m[10]; (*this)(2, 3) = m[14];
2108 (*this)(3, 0) = m[3]; (*this)(3, 1) = m[7]; (*this)(3, 2) = m[11]; (*this)(3, 3) = m[15];
2109 #endif
2110 }
2111
2112 /*----------------------------------------------------------------------------*/
2113 template <typename T>
2115 const Vec<4, T> &x,
2116 const Vec<4, T> &y,
2117 const Vec<4, T> &z,
2118 const Vec<4, T> &w
2119 ) {
2120
2121 #ifdef MATRIX_ROW_MAJOR
2122 (*this).set_row(0, x);
2123 (*this).set_row(1, y);
2124 (*this).set_row(2, z);
2125 (*this).set_row(3, w);
2126 #else
2127 (*this).set_col(0, x);
2128 (*this).set_col(1, y);
2129 (*this).set_col(2, z);
2130 (*this).set_col(3, w);
2131 #endif
2132 }
2133
2134 /*----------------------------------------------------------------------------*/
2135 template <typename T>
2136 inline Mat4<T>::Mat4(const Mat<3, 3, T> &rhs) {
2137 (*this)(0, 0) = rhs(0, 0); (*this)(0, 1) = rhs(0, 1); (*this)(0, 2) = rhs(0, 2); (*this)(0, 3) = T(0);
2138 (*this)(1, 0) = rhs(1, 0); (*this)(1, 1) = rhs(1, 1); (*this)(1, 2) = rhs(1, 2); (*this)(1, 3) = T(0);
2139 (*this)(2, 0) = rhs(2, 0); (*this)(2, 1) = rhs(2, 1); (*this)(2, 2) = rhs(2, 2); (*this)(2, 3) = T(0);
2140 (*this)(3, 0) = T(0); (*this)(3, 1) = T(0); (*this)(3, 2) = T(0); (*this)(3, 3) = T(1);
2141 }
2142
2143 /*----------------------------------------------------------------------------*/
2144 template <typename T>
2145 inline Mat3<T> Mat4<T>::sub() const {
2146 Mat3<T> mat;
2147 for (size_t i = 0; i < 3; i++) {
2148 for (size_t j = 0; j < 3; j++) {
2149 mat(i, j) = (*this)(i, j);
2150 }
2151 }
2152 return mat;
2153 }
2154
2155 /*----------------------------------------------------------------------------*/
2156 template <typename T>
2157 inline Mat4<T>::Mat4(const Vec<3, T> &s, const Quat<T> &rot, const Vec<3, T> &t) {
2158 assert(std::abs(rot.length() - 1) < epsilon<T>());
2159
2160 // get rotation matrix from quaternion:
2161 Mat3<T> r(rot);
2162
2163 // incorporate scale (cheaper than a matrix multiply):
2164 r(0, 0) = r(0, 0) * s[0]; r(0, 1) = r(0, 1) * s[1]; r(0, 2) = r(0, 2) * s[2];
2165 r(1, 0) = r(1, 0) * s[0]; r(1, 1) = r(1, 1) * s[1]; r(1, 2) = r(1, 2) * s[2];
2166 r(2, 0) = r(2, 0) * s[0]; r(2, 1) = r(2, 1) * s[1]; r(2, 2) = r(2, 2) * s[2];
2167
2168 // combine rotation matrix, scale and translation:
2169 (*this)(0, 0) = r(0, 0); (*this)(0, 1) = r(0, 1); (*this)(0, 2) = r(0, 2); (*this)(0, 3) = t[0];
2170 (*this)(1, 0) = r(1, 0); (*this)(1, 1) = r(1, 1); (*this)(1, 2) = r(1, 2); (*this)(1, 3) = t[1];
2171 (*this)(2, 0) = r(2, 0); (*this)(2, 1) = r(2, 1); (*this)(2, 2) = r(2, 2); (*this)(2, 3) = t[2];
2172 (*this)(3, 0) = T(0); (*this)(3, 1) = T(0); (*this)(3, 2) = T(0); (*this)(3, 3) = T(1);
2173 }
2174
2175 /*----------------------------------------------------------------------------*/
2176 template <typename T>
2178 return Mat4<T>(
2179 s, T(0), T(0), T(0),
2180 T(0), s, T(0), T(0),
2181 T(0), T(0), s, T(0),
2182 T(0), T(0), T(0), T(1)
2183 );
2184 }
2185
2186 /*----------------------------------------------------------------------------*/
2187 template <typename T>
2188 inline Mat4<T> Mat4<T>::scale(T x, T y, T z, T w) {
2189 return Mat4<T>(
2190 x, T(0), T(0), T(0),
2191 T(0), y, T(0), T(0),
2192 T(0), T(0), z, T(0),
2193 T(0), T(0), T(0), w
2194 );
2195 }
2196
2197 /*----------------------------------------------------------------------------*/
2198 template <typename T>
2200 return Mat4<T>(
2201 s.x, T(0), T(0), T(0),
2202 T(0), s.y, T(0), T(0),
2203 T(0), T(0), s.z, T(0),
2204 T(0), T(0), T(0), s.w
2205 );
2206 }
2207
2208 /*----------------------------------------------------------------------------*/
2209
2210 template <typename T>
2211 inline Mat4<T> Mat4<T>::rotation(const Vec<3, T> &axis, T angle) {
2212 assert(std::abs(axis.length() - 1) < epsilon<T>());
2213 return Mat4<T>(Mat3<T>::rotation(axis, angle)); // gen 3x3 rotation matrix as arg to Mat4 constructor
2214 }
2215
2216 /*----------------------------------------------------------------------------*/
2217
2218 template <typename T>
2219 inline Mat4<T> Mat4<T>::rotation(const Vec<3, T> &axis_angle) {
2220 const T len = axis_angle.length();
2221 return Mat4<T>(Mat3<T>::rotation(axis_angle/len, len));
2222 }
2223
2224 /*----------------------------------------------------------------------------*/
2225 template <typename T>
2227 // input must be unit quaternion
2228 assert(std::abs(q.length() - 1) < epsilon<T>());
2229 const T x = q.x;
2230 const T y = q.y;
2231 const T z = q.z;
2232 const T w = q.w;
2233 Mat4<T> m;
2234 m(0, 0) = 1.0 - 2.0 * (y*y + z*z); m(0, 1) = 2.0 * (x*y - w*z); m(0, 2) = 2.0 * (x*z + w*y); m(0, 3) = 0.0;
2235 m(1, 0) = 2.0 * (x*y + w*z); m(1, 1) = 1.0 - 2.0 * (x*x + z*z); m(1, 2) = 2.0 * (y*z - w*x); m(1, 3) = 0.0;
2236 m(2, 0) = 2.0 * (x*z - w*y); m(2, 1) = 2.0 * (y*z + w*x); m(2, 2) = 1.0 - 2.0 * (x*x + y*y); m(2, 3) = 0.0;
2237 m(3, 0) = 0.0; m(3, 1) = 0.0; m(3, 2) = 0.0; m(3, 3) = 1.0;
2238 return m;
2239 }
2240
2241 /*----------------------------------------------------------------------------*/
2242 template <typename T>
2243 inline Mat4<T> Mat4<T>::rotation(T x, T y, T z, int order) {
2244 // The code is not optimized. The final matrix can be directly derived.
2245 // http://www.songho.ca/opengl/gl_anglestoaxes.html
2246 // http://inside.mines.edu/fs_home/gmurray/ArbitraryAxisRotation/
2247 Mat4<T> rx; // Rotation about X-axis (Pitch)
2248 rx(0, 0) = T(1); rx(0, 1) = T(0); rx(0, 2) = T(0); rx(0, 3) = T(0);
2249 rx(1, 0) = T(0); rx(1, 1) = std::cos(x); rx(1, 2) = -std::sin(x); rx(1, 3) = T(0);
2250 rx(2, 0) = T(0); rx(2, 1) = std::sin(x); rx(2, 2) = std::cos(x); rx(2, 3) = T(0);
2251 rx(3, 0) = T(0); rx(3, 1) = T(0); rx(3, 2) = T(0); rx(3, 3) = T(1);
2252
2253 Mat4<T> ry; // Rotation about Y-axis (Yaw, Heading)
2254 ry(0, 0) = std::cos(y); ry(0, 1) = T(0); ry(0, 2) = std::sin(y); ry(0, 3) = T(0);
2255 ry(1, 0) = T(0); ry(1, 1) = T(1); ry(1, 2) = T(0); ry(1, 3) = T(0);
2256 ry(2, 0) = -std::sin(y); ry(2, 1) = T(0); ry(2, 2) = std::cos(y); ry(2, 3) = T(0);
2257 ry(3, 0) = T(0); ry(3, 1) = T(0); ry(3, 2) = T(0); ry(3, 3) = T(1);
2258
2259 Mat4<T> rz; // Rotation about Z-axis (Roll)
2260 rz(0, 0) = std::cos(z); rz(0, 1) = -std::sin(z); rz(0, 2) = T(0); rz(0, 3) = T(0);
2261 rz(1, 0) = std::sin(z); rz(1, 1) = std::cos(z); rz(1, 2) = T(0); rz(1, 3) = T(0);
2262 rz(2, 0) = T(0); rz(2, 1) = T(0); rz(2, 2) = T(1); rz(2, 3) = T(0);
2263 rz(3, 0) = T(0); rz(3, 1) = T(0); rz(3, 2) = T(0); rz(3, 3) = T(1);
2264
2265 switch (order) {
2266 case 123: return rz * ry * rx;
2267 case 132: return ry * rz * rx;
2268 case 213: return rz * rx * ry;
2269 case 231: return rx * rz * ry;
2270 case 312: return ry * rx * rz;
2271 case 321: return rx * ry * rz;
2272 default:
2273 LOG(ERROR) << "invalid rotation order";
2274 return rx * rz * ry;
2275 }
2276 }
2277
2278 /*----------------------------------------------------------------------------*/
2279 template <typename T>
2281 return Mat4<T>(
2282 T(1), T(0), T(0), t[0],
2283 T(0), T(1), T(0), t[1],
2284 T(0), T(0), T(1), t[2],
2285 T(0), T(0), T(0), T(1)
2286 );
2287 }
2288
2289 /*----------------------------------------------------------------------------*/
2290 template <typename T>
2291 inline Mat4<T> Mat4<T>::translation(T x, T y, T z) {
2292 return Mat4<T>(
2293 T(1), T(0), T(0), x,
2294 T(0), T(1), T(0), y,
2295 T(0), T(0), T(1), z,
2296 T(0), T(0), T(0), T(1)
2297 );
2298 }
2299
2300}
2301
2302#endif // EASY3D_CORE_MAT_H
2 by 2 matrix. Extends Mat with 2D-specific functionality and constructors.
Definition: mat.h:1460
static Mat2< T > scale(T s)
Static constructor return a 2D uniform scale matrix.
Definition: mat.h:1598
static Mat2< T > rotation(T angle)
Static constructor return a 2D rotation matrix.
Definition: mat.h:1589
Mat2()=default
Default constructor.
3x3 matrix. Extends Mat with 3D-specific functionality and constructors.
Definition: mat.h:1620
Mat3()=default
Default constructor.
static Mat3< T > rotation(const Vec< 3, T > &axis, T angle)
Static constructor returning a 3D rotation matrix defined by its axis and angle.
Definition: mat.h:1850
Mat2< T > sub() const
the upper-left 2x2 sub-matrix.
Definition: mat.h:1818
static Mat3< T > scale(T s)
Static constructor returning a 3D uniform scale matrix.
Definition: mat.h:1830
4x4 matrix. Extends Mat with 4D-specific functionality and constructors.
Definition: mat.h:1940
Mat4(const T *m)
Initialize elements from an array of type T.
Definition: mat.h:2096
Mat4()=default
Default constructor.
static Mat4< T > scale(const Vec< 4, T > &s)
Static constructor returning a 4D non-uniform scale matrix,.
Definition: mat.h:2199
Mat4(T s)
Initialized with diagonal as s and others zeros.
Definition: mat.h:2063
static Mat4< T > rotation(const Quat< T > &q)
Static constructor returning a 3D rotation matrix defined by a quaternion.
Definition: mat.h:2226
static Mat4< T > scale(T s)
Static constructor returning a 4D uniform scale matrix.
Definition: mat.h:2177
static Mat4< T > rotation(const Vec< 3, T > &axis, T angle)
Static constructor returning a 3D rotation matrix defined by its axis and angle.
Definition: mat.h:2211
Mat4(const Mat< 4, 4, T > &rhs)
Copy constructor. This provides compatibility with generic operations implemented by Mat.
Definition: mat.h:2075
Mat4(const Mat< 3, 3, T > &rhs)
Initialize from a 3D matrix. rhs forms the upper-left 3x3 sub-matrix, other elements are initialized ...
Definition: mat.h:2136
static Mat4< T > translation(const Vec< 3, T > &t)
Static constructor return a 3D translation matrix (as a 4D affine transformation).
Definition: mat.h:2280
Mat4(const Vec< 4, T > &x, const Vec< 4, T > &y, const Vec< 4, T > &z, const Vec< 4, T > &w)
Initialize elements from four vectors. If MATRIX_ROW_MAJOR is defined, x, y, z and w specify rows of ...
Definition: mat.h:2114
static Mat4< T > rotation(const Vec< 3, T > &axis_angle)
Static constructor returning a 3D rotation matrix defined by the axis–angle representation....
Definition: mat.h:2219
Mat4(const Vec< 3, T > &s, const Quat< T > &r, const Vec< 3, T > &t)
Initialize from scale/rotation/translation.
Definition: mat.h:2157
static Mat4< T > rotation(T x, T y, T z, int order=312)
Static constructor returning a 3D rotation matrix defined by Euler angles. The three rotations are ap...
Definition: mat.h:2243
Mat4(T s00, T s01, T s02, T s03, T s10, T s11, T s12, T s13, T s20, T s21, T s22, T s23, T s30, T s31, T s32, T s33)
Initialize elements from individual scalars. The digits following s in the parameter names indicate t...
Definition: mat.h:2082
Base class for matrix types.
Definition: mat.h:66
T & operator()(size_t r, size_t c)
Non-const row/column access to elements.
Definition: mat.h:613
Mat< N, M, T > & operator*=(T rhs)
Component-wise matrix-scalar multiplication/assignment.
Definition: mat.h:752
void swap_cols(size_t a, size_t b)
Swaps col a with col b.
Definition: mat.h:575
static Mat< N, M, T > identity()
Static constructor return an N x M identity matrix. see also load_identity()
Definition: mat.h:483
Mat()=default
Default constructor.
const T & operator()(size_t r, size_t c) const
Const row/column access to elements.
Definition: mat.h:600
Mat(const Mat< rN, rM, T > &rhs)
Copy constructor for rN >= N, rM >= M. For smaller incoming matrices (i.e, rN < N,...
Definition: mat.h:463
Mat(T s)
Initialized with diagonal as s and others zeros.
Definition: mat.h:450
Mat< N, M, T > & operator-=(const Mat< N, M, T > &rhs)
Component-wise matrix-matrix subtraction/assignment.
Definition: mat.h:741
Mat< N, M, T > & operator*=(const Mat< N, M, T > &rhs)
Matrix-matrix multiplication/assignment.
Definition: mat.h:724
Mat(const T *m)
Initialize elements from an array of type T.
Definition: mat.h:474
size_t num_rows() const
Return the number of rows.
Definition: mat.h:105
Mat< N, M, T > & operator-=(T rhs)
Component-wise matrix-scalar subtraction/assignment.
Definition: mat.h:776
Mat< N, M, T > operator*(T rhs) const
Component-wise matrix-scalar multiplication.
Definition: mat.h:704
void swap_rows(size_t a, size_t b)
Swaps row a with row b.
Definition: mat.h:542
void set_row(size_t r, const Vec< vN, T > &v)
Set row r from vector v. This copies the first M components from v, so vN must be >= M....
Definition: mat.h:521
Vec< M, T > row(size_t r) const
Return row r as a vector.
Definition: mat.h:498
void load_zero()
Set all elements 0.
Definition: mat.h:555
size_t num_columns() const
Return the number of columns.
Definition: mat.h:108
void load_identity(T s=T(1))
Set diagonal elements s and others 0.
Definition: mat.h:563
Mat< N, M, T > operator-(const Mat< N, M, T > &rhs) const
Component-wise matrix-matrix subtraction.
Definition: mat.h:667
Vec< N, T > col(size_t c) const
Return col c as a vector.
Definition: mat.h:509
bool operator!=(const Mat< N, M, T > &rhs) const
Inequality test.
Definition: mat.h:635
Mat< N, M, T > & operator/=(T rhs)
Component-wise matrix-scalar division/assignment.
Definition: mat.h:760
Mat< N, M, T > & operator+=(const Mat< N, M, T > &rhs)
Component-wise matrix-matrix addition/assignment.
Definition: mat.h:733
void set_col(size_t c, const Vec< vN, T > &v)
Set col c from vector v. This copies the first N components from v, so vN must be >= N....
Definition: mat.h:532
Vec< N, T > operator*(const Vec< M, T > &rhs) const
Matrix-vector multiplication. rhs must have the same number of elements as this matrix has columns....
Definition: mat.h:688
Mat< N, M, T > operator/(T rhs) const
Component-wise matrix-scalar division.
Definition: mat.h:713
Mat< N, M, T > & operator+=(T rhs)
Component-wise matrix-scalar addition/assignment.
Definition: mat.h:768
Mat< N, rM, T > operator*(const Mat< M, rM, T > &rhs) const
Matrix-matrix multiplication. rhs must have the same number of rows as this matrix has columns....
Definition: mat.h:643
bool operator==(const Mat< N, M, T > &rhs) const
Equality test.
Definition: mat.h:626
Mat< N, M, T > operator-() const
Component-wise matrix negation.
Definition: mat.h:676
Mat< N, M, T > operator+(const Mat< N, M, T > &rhs) const
Component-wise matrix-matrix addition.
Definition: mat.h:658
The Quaternion class represents 3D rotations and orientations.
Definition: quat.h:107
FT length() const
Return the length of the quaternion.
Definition: quat.h:285
A 2D vector (used for representing 2D points or vectors).
Definition: vec.h:320
A 3D vector (used for representing 3D points or vectors).
Definition: vec.h:439
A 4D vector (used for representing 3D points or vectors in homogeneous coordinates).
Definition: vec.h:584
Base class for vector types. It provides generic functionality for N dimensional vectors.
Definition: vec.h:34
T * data()
Returns the memory address of the vector.
Definition: vec.h:81
Definition: collider.cpp:182
Mat< M, N, T > transpose(const Mat< N, M, T > &m)
Transpose m.
Definition: mat.h:907
FT max()
Function returning maximum representable value for a given type.
T trace(const Mat< N, N, T > &m)
Return the trace (sum of elements on the main diagonal) of N by N (square) matrix m.
T determinant(const Mat< N, N, T > &m)
Return the determinant of N x N (square) matrix m.
bool has_nan(const GenericBox< DIM, FT > &box)
Definition: box.h:284
Mat< N, 1, FT > to_matrix(const Vec< N, FT > &v)
Convert a N-dimensional vector into a N by 1 matrix.
Definition: mat.h:1420
bool cholesky_decompose(const Mat< N, N, FT > &A, Mat< N, N, FT > &L)
Cholesky decomposition of a symmetric, positive definite matrix.
Definition: mat.h:1318
std::istream & operator>>(std::istream &is, GenericLine< DIM, FT > &line)
Definition: line.h:133
std::ostream & operator<<(std::ostream &os, Graph::Vertex v)
Definition: graph.h:920
void lu_back_substitution(const Mat< N, N, T > &alu, const Vec< N, T > &rowp, const Vec< N, T > &b, Vec< N, T > *x)
Solve a set of linear equations using outputs from lu_decomposition() as inputs.
Definition: mat.h:1277
T det(const Vec< 2, T > &v1, const Vec< 2, T > &v2)
Compute the determinant of the 2x2 matrix formed by the two 2D vectors as the two rows.
Definition: vec.h:405
bool lu_decomposition(const Mat< N, N, T > &a, Mat< N, N, T > *alu, Vec< N, T > *rowp, T *d)
Perform LU decomposition of a square matrix.
Definition: mat.h:1201
Mat< N, N, T > inverse(const Mat< N, N, T > &m)
Return the inverse of N x N (square) matrix m.
Definition: mat.h:977
std::vector< FT > sum(const Matrix< FT > &)
Definition: matrix.h:1454
bool gauss_jordan_elimination(const Mat< N, N, T > &a, const Mat< N, M, T > &b, Mat< N, N, T > *ainv, Mat< N, M, T > *x)
Perform Gauss-Jordan elimination to solve a set of linear equations and additionally compute the inve...
Definition: mat.h:1121
void cholesky_solve(const Mat< N, N, FT > &L, const Vec< N, FT > &b, Vec< N, FT > &x)
Definition: mat.h:1345
Mat< N, M, T > operator*(T lhs, const Mat< N, M, T > &rhs)
Component-wise scalar-matrix multiplication.
Definition: mat.h:787
Mat< N, M, T > tensor(const Vec< M, T > &u, const Vec< N, T > &v)
Definition: mat.h:1112