MorphoGraphX  2.0-1-227
Matrix.hpp
Go to the documentation of this file.
1 //
2 // This file is part of MorphoGraphX - http://www.MorphoGraphX.org
3 // Copyright (C) 2012-2015 Richard S. Smith and collaborators.
4 //
5 // If you use MorphoGraphX in your work, please cite:
6 // http://dx.doi.org/10.7554/eLife.05864
7 //
8 // MorphoGraphX is free software, and is licensed under under the terms of the
9 // GNU General (GPL) Public License version 2.0, http://www.gnu.org/licenses.
10 //
11 #ifndef MATRIX_HPP
12 #define MATRIX_HPP
13 
20 #include <Config.hpp>
21 
22 #include <cuda/CudaGlobal.hpp>
23 #include <StaticAssert.hpp>
24 #include <Util.hpp>
25 #include <Vector.hpp>
26 
27 #include <cmath>
28 #ifndef __CUDACC__
29 # include <QTextStream>
30 #endif
31 
32 namespace mgx
33 {
39  template <size_t nRows, size_t nCols, typename T = double> class Matrix {
40  private:
41  Vector<nCols, T> rows[nRows];
42 
43  public:
44  typedef T value_type;
45  typedef T& reference_type;
46  typedef const T& const_reference_type;
47  typedef T* pointer_type;
48  typedef const T* const_pointer_type;
49  typedef T* iterator;
50  typedef const T* const_iterator;
51 
52  static const size_t numcols = nCols;
53  static const size_t numrows = nRows;
54 
59  Matrix(void)
60  {
61  for(size_t i = 0; i < nRows; i++)
62  for(size_t j = 0; j < nCols; j++)
63  rows[i][j] = 0;
64  }
65 
71  template <typename T1> CU_HOST_DEVICE Matrix(const Matrix<nRows, nCols, T1>& mat)
72  {
73  for(size_t i = 0; i < nRows; i++)
74  for(size_t j = 0; j < nCols; j++)
75  rows[i][j] = T(mat(i, j));
76  }
77 
83  template <typename T1> CU_HOST_DEVICE Matrix(const Vector<nCols, T1>* vecs)
84  {
85  for(size_t i = 0; i < nRows; i++)
86  rows[i] = vecs[i];
87  }
88 
92  template <typename T1> CU_HOST_DEVICE
93  Matrix(const Vector<nCols,T1> &r0, const Vector<nCols,T1> &r1, const Vector<nCols,T1> &r2)
94  {
95  STATIC_ASSERT(nRows == 3);
96  rows[0] = r0;
97  rows[1] = r1;
98  rows[2] = r2;
99  }
100 
110  template <typename T1> CU_HOST_DEVICE
111  Matrix(const T1* values)
112  {
113  for(size_t i = 0; i < nRows; i++)
114  rows[i] = Vector<nCols, T>(values + (i * nCols));
115  }
116 
127  Matrix(const T* values)
128  {
129  for(size_t i = 0; i < nRows; i++)
130  rows[i] = Vector<nCols, T>(values + (i * nCols));
131  }
132 
139  Matrix(const T& value)
140  {
141  for(size_t i = 0; i < nRows; i++)
142  for(size_t j = 0; j < nCols; j++)
143  rows[i][j] = (i == j) ? value : 0;
144  }
145 
153  return Vector<2, size_t>(nRows, nCols);
154  }
155 
160  static size_t nbRows() {
161  return nRows;
162  }
163 
168  static size_t nbColumns() {
169  return nCols;
170  }
171 
178  const T* c_data() const {
179  return rows[0].c_data();
180  }
181 
188  T* data() {
189  return rows[0].data();
190  }
191 
196  Matrix operator-(void) const
197  {
198  Matrix ans;
199 
200  for(size_t i = 0; i < nRows; i++)
201  ans.rows[i] = -rows[i];
202 
203  return ans;
204  }
205 
210  Matrix operator+(const Matrix& mat) const
211  {
212  Matrix ans;
213 
214  for(size_t i = 0; i < nRows; i++)
215  ans.rows[i] = rows[i] + mat.rows[i];
216 
220  return ans;
221  }
222 
224  Matrix operator-(const Matrix& mat) const
225  {
226  Matrix ans;
227 
228  for(size_t i = 0; i < nRows; i++)
229  ans.rows[i] = rows[i] - mat.rows[i];
230 
231  return ans;
232  }
233 
238  Matrix operator*(const T& scalar) const
239  {
240  Matrix ans;
241 
242  for(size_t i = 0; i < nRows; i++)
243  ans[i] = rows[i] * scalar;
244 
245  return ans;
246  }
247 
252  Matrix operator/(const T& scalar) const
253  {
254  Matrix ans;
255 
256  for(size_t i = 0; i < nRows; i++)
257  ans[i] = rows[i] / scalar;
258 
259  return ans;
260  }
261 
266  friend Matrix operator*(const T& scalar, const Matrix& mat)
267  {
268  Matrix ans;
269 
270  for(size_t i = 0; i < nRows; i++)
271  ans[i] = scalar * mat.rows[i];
272 
273  return ans;
274  }
275 
281  {
282  Vector<nRows, T> ans;
283  for(size_t i = 0; i < nRows; ++i) {
284  T value = 0;
285  for(size_t j = 0; j < nCols; ++j)
286  value += rows[i][j] * vec[j];
287  ans[i] = value;
288  }
289  return ans;
290  }
291 
293  Matrix& operator=(const Matrix& mat)
294  {
295  for(size_t i = 0; i < nRows; i++)
296  rows[i] = mat.rows[i];
297 
298  return (*this);
299  }
300 
302  Matrix& operator+=(const Matrix& mat) {
303  return ((*this) = (*this) + mat);
304  }
305 
307  Matrix& operator-=(const Matrix& mat) {
308  return ((*this) = (*this) - mat);
309  }
310 
312  Matrix& operator*=(const T& scalar) {
313  return ((*this) = (*this) * scalar);
314  }
315 
317  Matrix& operator/=(const T& scalar) {
318  return ((*this) = (*this) / scalar);
319  }
320 
322  Matrix& operator*=(const Matrix& mat) {
323  return ((*this) = (*this) * mat);
324  }
325 
330  bool operator==(const Matrix& mat) const
331  {
332  for(size_t i = 0; i < nRows; i++)
333  if(rows[i] != mat.rows[i])
334  return false;
335 
336  return true;
337  }
338 
343  bool operator!=(const Matrix& mat) const {
344  return (!((*this) == mat));
345  }
346 
353  bool operator<(const Matrix& other) const
354  {
355  for(size_t i = 0; i < nRows; ++i) {
356  if(rows[i] < other.rows[i])
357  return true;
358  if(rows[i] > other.rows[i])
359  return false;
360  }
361  return false;
362  }
363 
364  #ifndef __CUDACC__
365  friend QTextStream& operator<<(QTextStream& out, const Matrix& mat)
366  {
367  for(size_t i = 0; i < nRows; i++) {
368  out << mat.rows[i];
369  if(i != (nRows - 1))
370  out << " ";
371  }
372 
373  return out;
374  }
375 
376  friend QTextStream& operator>>(QTextStream& in, Matrix& mat)
377  {
378  for(size_t i = 0; i < nRows && !in.atEnd(); i++)
379  in >> mat.rows[i];
380  return in;
381  }
382  #endif
383 
385  friend std::ostream& operator<<(std::ostream& out, const Matrix& mat)
386  {
387  for(size_t i = 0; i < nRows; i++) {
388  out << mat.rows[i];
389  if(i != (nRows - 1))
390  out << " ";
391  }
392 
393  return out;
394  }
395 
397  friend std::istream& operator>>(std::istream& in, Matrix& mat)
398  {
399  for(size_t i = 0; i < nRows && in; i++)
400  in >> mat.rows[i];
401  return in;
402  }
403 
411  return rows[idx];
412  }
413 
420  Vector<nCols, T> operator[](size_t idx) const {
421  return rows[idx];
422  }
423 
428  T& operator()(size_t i, size_t j) {
429  return rows[i][j];
430  }
431 
436  T operator()(size_t i, size_t j) const {
437  return rows[i][j];
438  }
439 
444  static Matrix identity()
445  {
446  Matrix mat(1);
447  return mat;
448  }
449 
454  Matrix& zero(void)
455  {
456  for(size_t i = 0; i < nRows; i++)
457  for(size_t j = 0; j < nCols; j++)
458  rows[i][j] = 0.0;
459  return (*this);
460  }
461 
468  Matrix& operator=(const T& value)
469  {
470  for(size_t i = 0; i < nRows; ++i) {
471  for(size_t j = 0; j < nCols; ++j) {
472  if(i == j)
473  rows[i][j] = value;
474  else
475  rows[i][j] = 0;
476  }
477  }
478  return *this;
479  }
480 
486  {
488  for(size_t i = 0; i < nRows; ++i)
489  for(size_t j = 0; j < nCols; ++j)
490  t[i][j] = rows[j][i];
491  return t;
492  }
493 
501  static Matrix<3, 3, T> rotation(const Vector<3, T>& direction, T angle)
502  {
503  T ca = std::cos(angle);
504  T sa = std::sin(angle);
505  Matrix<3, 3, T> r;
506  T x = direction.x();
507  T y = direction.y();
508  T z = direction.z();
509  r[0].set(ca + (1 - ca) * x * x, (1 - ca) * x * y - sa * z, (1 - ca) * z * x + sa * y);
510  r[1].set((1 - ca) * y * x + sa * z, ca + (1 - ca) * y * y, (1 - ca) * z * y - sa * x);
511  r[2].set((1 - ca) * x * z - sa * y, (1 - ca) * y * z + sa * x, ca + (1 - ca) * z * z);
512  return r;
513  }
514 
522  static Matrix<4, 4, T> rotation(const Vector<4, T>& direction, T angle)
523  {
524  T ca = std::cos(angle);
525  T sa = std::sin(angle);
526  Matrix<4, 4, T> r;
527  T x = direction.x() / direction.t();
528  T y = direction.y() / direction.t();
529  T z = direction.z() / direction.t();
530  r[0].set(ca + (1 - ca) * x * x, (1 - ca) * x * y - sa * z, (1 - ca) * z * x + sa * y, 0);
531  r[1].set((1 - ca) * y * x + sa * z, ca + (1 - ca) * y * y, (1 - ca) * z * y - sa * x, 0);
532  r[2].set((1 - ca) * x * z - sa * y, (1 - ca) * y * z + sa * x, ca + (1 - ca) * z * z, 0);
533  r[3].set(0, 0, 0, 1);
534  return r;
535  }
536 
541  T trace() const
542  {
543  T acc = 0;
544  for(size_t i = 0; i < min(nRows, nCols); ++i) {
545  acc += rows[i][i];
546  }
547  return acc;
548  }
549 
551  void fillArray(T* array)
552  {
553  memcpy(array, &rows[0][0], sizeof(T) * nRows * nCols);
554  }
555 
561  {
562  STATIC_ASSERT(nRows == nCols);
563  Vector<nRows, T> res;
564  for(size_t i = 0; i < nRows; ++i)
565  res[i] = (*this)(i, i);
566  return res;
567  }
568  };
569 
574  template <size_t nRows, size_t nCols, typename T>
575  CU_HOST_DEVICE Vector<nCols, T> operator*(const Vector<nCols, T>& vec, const Matrix<nRows, nCols, T>& mat);
576 
581  template <size_t nRows, size_t nSize, size_t nCols, typename T>
582  CU_HOST_DEVICE Matrix<nRows, nCols, T> operator*(const Matrix<nRows, nSize, T>& mat1,
583  const Matrix<nSize, nCols, T>& mat2);
584 
589  template <typename T> T CU_HOST_DEVICE det(const Matrix<1, 1, T>& mat);
594  template <typename T> T CU_HOST_DEVICE det(const Matrix<2, 2, T>& mat);
599  template <typename T> T CU_HOST_DEVICE det(const Matrix<3, 3, T>& mat);
606  template <size_t nRows, typename T> T CU_HOST_DEVICE det(const Matrix<nRows, nRows, T>& mat);
607 
608  template <size_t nRows, typename T>
609  CU_HOST_DEVICE T matrix_minor(const Matrix<nRows, nRows, T>& mat, size_t i, size_t j);
610 
615  template <size_t nRows, typename T> CU_HOST_DEVICE T cofactor(const Matrix<nRows, nRows, T>& mat, size_t i, size_t j);
616 
622  template <typename T> CU_HOST_DEVICE Matrix<1, 1, T> inverse(const Matrix<1, 1, T>& mat);
628  template <typename T> CU_HOST_DEVICE Matrix<2, 2, T> inverse(const Matrix<2, 2, T>& mat);
634  template <typename T> CU_HOST_DEVICE Matrix<3, 3, T> inverse(const Matrix<3, 3, T>& mat);
640  template <typename T> CU_HOST_DEVICE Matrix<4, 4, T> inverse(const Matrix<4, 4, T>& mat);
641 
649  template <size_t nRows, typename T> CU_HOST_DEVICE Matrix<nRows, nRows, T> inverse(const Matrix<nRows, nRows, T>& mat);
650 
656  template <size_t nRows, size_t nCols, typename T>
657  CU_HOST_DEVICE Matrix<nCols, nRows, T> transpose(const Matrix<nRows, nCols, T>& mat);
658 
667  template <size_t nRows, size_t nCols, typename T> CU_HOST_DEVICE T norm(const Matrix<nRows, nCols, T>& mat);
668 
676  template <size_t nRows, size_t nCols, typename T> CU_HOST_DEVICE T normsq(const Matrix<nRows, nCols, T>& mat);
677 
683  template <size_t nRows, size_t nCols, typename T>
685  {
687  for(size_t i = 0; i < nRows; ++i) {
688  const Vector<nCols, T>& mrow = m[i];
689  Vector<nCols, T>& rrow = result[i];
690  for(size_t j = 0; j < nCols; ++j) {
691  rrow[j] = (*fct)(mrow[j]);
692  }
693  }
694  return result;
695  }
696 
702  template <size_t nRows, size_t nCols, typename T>
704  {
706  for(size_t i = 0; i < nRows; ++i) {
707  const Vector<nCols, T>& mrow = m[i];
708  Vector<nCols, T>& rrow = result[i];
709  for(size_t j = 0; j < nCols; ++j) {
710  rrow[j] = (*fct)(mrow[j]);
711  }
712  }
713  return result;
714  }
715 
721  template <size_t nRows, size_t nCols, typename T>
723  {
725  for(size_t i = 0; i < nRows; ++i) {
726  const Vector<nCols, T>& mrow = m[i];
727  Vector<nCols, T>& rrow = result[i];
728  for(size_t j = 0; j < nCols; ++j) {
729  rrow[j] = (*fct)(mrow[j]);
730  }
731  }
732  return result;
733  }
734 
740  template <size_t nRows, size_t nCols, typename T, typename T1>
742  {
744  for(size_t i = 0; i < nRows; ++i) {
745  const Vector<nCols, T1>& mrow = m[i];
746  Vector<nCols, T>& rrow = result[i];
747  for(size_t j = 0; j < nCols; ++j) {
748  rrow[j] = (*fct)(mrow[j]);
749  }
750  }
751  return result;
752  }
753 
759  template <size_t nRows, size_t nCols, typename T, typename T1>
761  {
763  for(size_t i = 0; i < nRows; ++i) {
764  const Vector<nCols, T1>& mrow = m[i];
765  Vector<nCols, T>& rrow = result[i];
766  for(size_t j = 0; j < nCols; ++j) {
767  rrow[j] = (*fct)(mrow[j]);
768  }
769  }
770  return result;
771  }
772 
778  template <size_t nRows, size_t nCols, typename T, typename T1>
780  {
782  for(size_t i = 0; i < nRows; ++i) {
783  const Vector<nCols, T1>& mrow = m[i];
784  Vector<nCols, T>& rrow = result[i];
785  for(size_t j = 0; j < nCols; ++j) {
786  rrow[j] = (*fct)(mrow[j]);
787  }
788  }
789  return result;
790  }
791 
797  template <size_t nRows, size_t nCols, typename T>
799  const Matrix<nRows, nCols, T>& m2)
800  {
802  for(size_t i = 0; i < nRows; ++i) {
803  const Vector<nCols, T>& mrow1 = m1[i];
804  const Vector<nCols, T>& mrow2 = m2[i];
805  Vector<nCols, T>& rrow = result[i];
806  for(size_t j = 0; j < nCols; ++j) {
807  rrow[j] = (*fct)(mrow1[j], mrow2[j]);
808  }
809  }
810  return result;
811  }
812 
818  template <size_t nRows, size_t nCols, typename T>
819  CU_HOST_DEVICE Matrix<nRows, nCols, T> map(T (*fct)(const T&, const T&), const Matrix<nRows, nCols, T>& m1,
820  const Matrix<nRows, nCols, T>& m2)
821  {
823  for(size_t i = 0; i < nRows; ++i) {
824  const Vector<nCols, T>& mrow1 = m1[i];
825  const Vector<nCols, T>& mrow2 = m2[i];
826  Vector<nCols, T>& rrow = result[i];
827  for(size_t j = 0; j < nCols; ++j) {
828  rrow[j] = (*fct)(mrow1[j], mrow2[j]);
829  }
830  }
831  return result;
832  }
833 
839  template <size_t nRows, size_t nCols, typename T>
840  CU_HOST_DEVICE Matrix<nRows, nCols, T> map(const T& (*fct)(const T &, const T &), const Matrix<nRows, nCols, T>& m1,
841  const Matrix<nRows, nCols, T>& m2)
842  {
844  for(size_t i = 0; i < nRows; ++i) {
845  const Vector<nCols, T>& mrow1 = m1[i];
846  const Vector<nCols, T>& mrow2 = m2[i];
847  Vector<nCols, T>& rrow = result[i];
848  for(size_t j = 0; j < nCols; ++j) {
849  rrow[j] = (*fct)(mrow1[j], mrow2[j]);
850  }
851  }
852  return result;
853  }
854 
860  template <size_t nRows, size_t nCols, typename T, typename T1, typename T2>
862  const Matrix<nRows, nCols, T2>& m2)
863  {
865  for(size_t i = 0; i < nRows; ++i) {
866  const Vector<nCols, T1>& mrow1 = m1[i];
867  const Vector<nCols, T2>& mrow2 = m2[i];
868  Vector<nCols, T>& rrow = result[i];
869  for(size_t j = 0; j < nCols; ++j) {
870  rrow[j] = (*fct)(mrow1[j], mrow2[j]);
871  }
872  }
873  return result;
874  }
875 
881  template <size_t nRows, size_t nCols, typename T, typename T1, typename T2>
882  CU_HOST_DEVICE Matrix<nRows, nCols, T> map(T (*fct)(const T1&, const T2&), const Matrix<nRows, nCols, T1>& m1,
883  const Matrix<nRows, nCols, T2>& m2)
884  {
886  for(size_t i = 0; i < nRows; ++i) {
887  const Vector<nCols, T1>& mrow1 = m1[i];
888  const Vector<nCols, T2>& mrow2 = m2[i];
889  Vector<nCols, T>& rrow = result[i];
890  for(size_t j = 0; j < nCols; ++j) {
891  rrow[j] = (*fct)(mrow1[j], mrow2[j]);
892  }
893  }
894  return result;
895  }
896 
902  template <size_t nRows, size_t nCols, typename T, typename T1, typename T2>
903  CU_HOST_DEVICE Matrix<nRows, nCols, T> map(const T& (*fct)(const T1 &, const T2 &), const Matrix<nRows, nCols, T1>& m1,
904  const Matrix<nRows, nCols, T2>& m2)
905  {
907  for(size_t i = 0; i < nRows; ++i) {
908  const Vector<nCols, T1>& mrow1 = m1[i];
909  const Vector<nCols, T2>& mrow2 = m2[i];
910  Vector<nCols, T>& rrow = result[i];
911  for(size_t j = 0; j < nCols; ++j) {
912  rrow[j] = (*fct)(mrow1[j], mrow2[j]);
913  }
914  }
915  return result;
916  }
917 
918  template <size_t nRows, size_t nCols, typename T>
920  {
921  Vector<nCols, T> ans;
922  for(size_t i = 0; i < nCols; ++i) {
923  T value = 0;
924  for(size_t j = 0; j < nRows; ++j)
925  value += mat(j, i) * vec[j];
926  ans[i] = value;
927  }
928  return ans;
929  }
930 
931  template <size_t nRows, size_t intSize, size_t nCols, typename T>
934  {
936  for(size_t i = 0; i < nRows; i++)
937  for(size_t j = 0; j < nCols; j++) {
938  T acc = 0;
939  for(size_t k = 0; k < intSize; k++)
940  acc += mat1(i, k) * mat2(k, j);
941  ans(i, j) = acc;
942  }
943  return ans;
944  }
945 
946  template <typename T> CU_HOST_DEVICE T det(const Matrix<1, 1, T>& mat) {
947  return mat(0, 0);
948  }
949 
950  template <typename T> CU_HOST_DEVICE T det(const Matrix<2, 2, T>& mat)
951  {
952  return mat(0, 0) * mat(1, 1) - mat(0, 1) * mat(1, 0);
953  }
954 
955  template <typename T> CU_HOST_DEVICE T det(const Matrix<3, 3, T>& mat)
956  {
957  return mat(0, 0) * mat(1, 1) * mat(2, 2) + mat(0, 1) * mat(1, 2) * mat(2, 0) + mat(0, 2) * mat(1, 0) * mat(2, 1)
958  - mat(0, 0) * mat(1, 2) * mat(2, 1) - mat(0, 1) * mat(1, 0) * mat(2, 2) - mat(0, 2) * mat(1, 1) * mat(2, 0);
959  }
960 
961  template <size_t nRows, typename T> CU_HOST_DEVICE T det(const Matrix<nRows, nRows, T>& mat)
962  {
963  STATIC_ASSERT(nRows > 3);
964  T acc = 0;
965  for(size_t i = 0; i < nRows; i++) {
966  acc += mat(i, 0) * cofactor(mat, i, 0);
967  }
968  return acc;
969  }
970 
971  template <size_t nRows, typename T>
972  CU_HOST_DEVICE T matrix_minor(const Matrix<nRows, nRows, T>& mat, size_t i, size_t j)
973  {
974  STATIC_ASSERT(nRows > 0);
975  Matrix<nRows - 1, nRows - 1, T> ans;
976  for(size_t i1 = 0, i2 = 0; i1 < nRows; i1++, i2++) {
977  if(i1 == i) {
978  i2--;
979  } else {
980  for(size_t j1 = 0, j2 = 0; j1 < nRows; j1++, j2++) {
981  if(j1 == j) {
982  j2--;
983  } else {
984  ans(i2, j2) = mat(i1, j1);
985  }
986  }
987  }
988  }
989  return det(ans);
990  }
991 
992  template <size_t nRows, typename T> CU_HOST_DEVICE T cofactor(const Matrix<nRows, nRows, T>& mat, size_t i, size_t j)
993  {
994  T inv = 1;
995  if((i + j) % 2) {
996  inv = -1;
997  }
998  return inv * matrix_minor(mat, i, j);
999  }
1000 
1001  template <typename T> CU_HOST_DEVICE Matrix<1, 1, T> inverse(const Matrix<1, 1, T>& mat)
1002  {
1003  Matrix<1, 1, T> ans;
1004  ans[0][0] = 1 / mat(0, 0);
1005  return ans;
1006  }
1007 
1008  template <typename T> CU_HOST_DEVICE Matrix<2, 2, T> inverse(const Matrix<2, 2, T>& mat)
1009  {
1010  Matrix<2, 2, T> ans;
1011  T d;
1012  d = det(mat);
1013  if(d == 0)
1014  return ans;
1015  ans(0, 0) = mat(1, 1) / d;
1016  ans(0, 1) = mat(0, 1) / -d;
1017  ans(1, 0) = mat(1, 0) / -d;
1018  ans(1, 1) = mat(0, 0) / d;
1019  return ans;
1020  }
1021 
1022  template <typename T> CU_HOST_DEVICE Matrix<3, 3, T> inverse(const Matrix<3, 3, T>& mat)
1023  {
1024  Matrix<3, 3, T> ans;
1025  T d;
1026  d = det(mat);
1027  if(d == 0)
1028  return ans;
1029  ans(0, 0) = (mat(1, 1) * mat(2, 2) - mat(1, 2) * mat(2, 1)) / d;
1030  ans(1, 1) = (mat(0, 0) * mat(2, 2) - mat(0, 2) * mat(2, 0)) / d;
1031  ans(2, 2) = (mat(1, 1) * mat(0, 0) - mat(1, 0) * mat(0, 1)) / d;
1032  ans(1, 0) = (mat(1, 2) * mat(2, 0) - mat(1, 0) * mat(2, 2)) / d;
1033  ans(0, 1) = (mat(2, 1) * mat(0, 2) - mat(0, 1) * mat(2, 2)) / d;
1034  ans(2, 0) = (mat(1, 0) * mat(2, 1) - mat(1, 1) * mat(2, 0)) / d;
1035  ans(0, 2) = (mat(0, 1) * mat(1, 2) - mat(1, 1) * mat(0, 2)) / d;
1036  ans(1, 2) = (mat(0, 2) * mat(1, 0) - mat(0, 0) * mat(1, 2)) / d;
1037  ans(2, 1) = (mat(2, 0) * mat(0, 1) - mat(0, 0) * mat(2, 1)) / d;
1038  return ans;
1039  }
1040 
1041  template <typename T> Matrix<4, 4, T> inverse(const Matrix<4, 4, T>& m)
1042  {
1043  Matrix<4, 4, T> inv;
1044 
1045  inv(0, 0) = (m(1, 1) * m(2, 2) * m(3, 3) - m(1, 1) * m(3, 2) * m(2, 3) - m(1, 2) * m(2, 1) * m(3, 3)
1046  + m(1, 2) * m(3, 1) * m(2, 3) + m(1, 3) * m(2, 1) * m(3, 2) - m(1, 3) * m(3, 1) * m(2, 2));
1047 
1048  inv(0, 1) = (-m(0, 1) * m(2, 2) * m(3, 3) + m(0, 1) * m(3, 2) * m(2, 3) + m(0, 2) * m(2, 1) * m(3, 3)
1049  - m(0, 2) * m(3, 1) * m(2, 3) - m(0, 3) * m(2, 1) * m(3, 2) + m(0, 3) * m(3, 1) * m(2, 2));
1050 
1051  inv(0, 2) = (m(0, 1) * m(1, 2) * m(3, 3) - m(0, 1) * m(3, 2) * m(1, 3) - m(0, 2) * m(1, 1) * m(3, 3)
1052  + m(0, 2) * m(3, 1) * m(1, 3) + m(0, 3) * m(1, 1) * m(3, 2) - m(0, 3) * m(3, 1) * m(1, 2));
1053 
1054  inv(0, 3) = (-m(0, 1) * m(1, 2) * m(2, 3) + m(0, 1) * m(2, 2) * m(1, 3) + m(0, 2) * m(1, 1) * m(2, 3)
1055  - m(0, 2) * m(2, 1) * m(1, 3) - m(0, 3) * m(1, 1) * m(2, 2) + m(0, 3) * m(2, 1) * m(1, 2));
1056 
1057  inv(1, 0) = (-m(1, 0) * m(2, 2) * m(3, 3) + m(1, 0) * m(3, 2) * m(2, 3) + m(1, 2) * m(2, 0) * m(3, 3)
1058  - m(1, 2) * m(3, 0) * m(2, 3) - m(1, 3) * m(2, 0) * m(3, 2) + m(1, 3) * m(3, 0) * m(2, 2));
1059 
1060  inv(1, 1) = (m(0, 0) * m(2, 2) * m(3, 3) - m(0, 0) * m(3, 2) * m(2, 3) - m(0, 2) * m(2, 0) * m(3, 3)
1061  + m(0, 2) * m(3, 0) * m(2, 3) + m(0, 3) * m(2, 0) * m(3, 2) - m(0, 3) * m(3, 0) * m(2, 2));
1062 
1063  inv(1, 2) = (-m(0, 0) * m(1, 2) * m(3, 3) + m(0, 0) * m(3, 2) * m(1, 3) + m(0, 2) * m(1, 0) * m(3, 3)
1064  - m(0, 2) * m(3, 0) * m(1, 3) - m(0, 3) * m(1, 0) * m(3, 2) + m(0, 3) * m(3, 0) * m(1, 2));
1065 
1066  inv(1, 3) = (m(0, 0) * m(1, 2) * m(2, 3) - m(0, 0) * m(2, 2) * m(1, 3) - m(0, 2) * m(1, 0) * m(2, 3)
1067  + m(0, 2) * m(2, 0) * m(1, 3) + m(0, 3) * m(1, 0) * m(2, 2) - m(0, 3) * m(2, 0) * m(1, 2));
1068 
1069  inv(2, 0) = (m(1, 0) * m(2, 1) * m(3, 3) - m(1, 0) * m(3, 1) * m(2, 3) - m(1, 1) * m(2, 0) * m(3, 3)
1070  + m(1, 1) * m(3, 0) * m(2, 3) + m(1, 3) * m(2, 0) * m(3, 1) - m(1, 3) * m(3, 0) * m(2, 1));
1071 
1072  inv(2, 1) = (-m(0, 0) * m(2, 1) * m(3, 3) + m(0, 0) * m(3, 1) * m(2, 3) + m(0, 1) * m(2, 0) * m(3, 3)
1073  - m(0, 1) * m(3, 0) * m(2, 3) - m(0, 3) * m(2, 0) * m(3, 1) + m(0, 3) * m(3, 0) * m(2, 1));
1074 
1075  inv(2, 2) = (m(0, 0) * m(1, 1) * m(3, 3) - m(0, 0) * m(3, 1) * m(1, 3) - m(0, 1) * m(1, 0) * m(3, 3)
1076  + m(0, 1) * m(3, 0) * m(1, 3) + m(0, 3) * m(1, 0) * m(3, 1) - m(0, 3) * m(3, 0) * m(1, 1));
1077 
1078  inv(2, 3) = (-m(0, 0) * m(1, 1) * m(2, 3) + m(0, 0) * m(2, 1) * m(1, 3) + m(0, 1) * m(1, 0) * m(2, 3)
1079  - m(0, 1) * m(2, 0) * m(1, 3) - m(0, 3) * m(1, 0) * m(2, 1) + m(0, 3) * m(2, 0) * m(1, 1));
1080 
1081  inv(3, 0) = (-m(1, 0) * m(2, 1) * m(3, 2) + m(1, 0) * m(3, 1) * m(2, 2) + m(1, 1) * m(2, 0) * m(3, 2)
1082  - m(1, 1) * m(3, 0) * m(2, 2) - m(1, 2) * m(2, 0) * m(3, 1) + m(1, 2) * m(3, 0) * m(2, 1));
1083 
1084  inv(3, 1) = (m(0, 0) * m(2, 1) * m(3, 2) - m(0, 0) * m(3, 1) * m(2, 2) - m(0, 1) * m(2, 0) * m(3, 2)
1085  + m(0, 1) * m(3, 0) * m(2, 2) + m(0, 2) * m(2, 0) * m(3, 1) - m(0, 2) * m(3, 0) * m(2, 1));
1086 
1087  inv(3, 2) = (-m(0, 0) * m(1, 1) * m(3, 2) + m(0, 0) * m(3, 1) * m(1, 2) + m(0, 1) * m(1, 0) * m(3, 2)
1088  - m(0, 1) * m(3, 0) * m(1, 2) - m(0, 2) * m(1, 0) * m(3, 1) + m(0, 2) * m(3, 0) * m(1, 1));
1089 
1090  inv(3, 3) = (m(0, 0) * m(1, 1) * m(2, 2) - m(0, 0) * m(2, 1) * m(1, 2) - m(0, 1) * m(1, 0) * m(2, 2)
1091  + m(0, 1) * m(2, 0) * m(1, 2) + m(0, 2) * m(1, 0) * m(2, 1) - m(0, 2) * m(2, 0) * m(1, 1));
1092 
1093  T det = m(0, 0) * inv(0, 0) + m(1, 0) * inv(0, 1) + m(2, 0) * inv(0, 2) + m(3, 0) * inv(0, 3);
1094 
1095  T inv_det = 1.0 / det;
1096 
1097  return inv_det * inv;
1098  }
1099 
1100  template <size_t nRows, typename T> CU_HOST_DEVICE Matrix<nRows, nRows, T> inverse(const Matrix<nRows, nRows, T>& mat)
1101  {
1103  T d = det(mat);
1104  if(d == 0)
1105  return ans;
1106  for(size_t i = 0; i < nRows; ++i)
1107  for(size_t j = 0; j < nRows; ++j) {
1108  ans(i, j) = cofactor(mat, j, i) / d;
1109  }
1110  return ans;
1111  }
1112 
1113  template <size_t nRows, size_t nCols, typename T>
1115  {
1117  for(size_t i = 0; i < nRows; i++)
1118  for(size_t j = 0; j < nCols; j++)
1119  ans(j, i) = mat(i, j);
1120  return ans;
1121  }
1122 
1123  template <size_t nRows, size_t nCols, typename T> CU_HOST_DEVICE T normsq(const Matrix<nRows, nCols, T>& mat)
1124  {
1125  T acc = 0;
1126  for(size_t i = 0; i < nRows; i++)
1127  for(size_t j = 0; j < nCols; j++) {
1128  const T& v = mat(i, j);
1129  acc += v * v;
1130  }
1131  return acc;
1132  }
1133 
1134  template <size_t nRows, size_t nCols, typename T> CU_HOST_DEVICE T norm(const Matrix<nRows, nCols, T>& mat)
1135  {
1136  return std::sqrt(normsq(mat));
1137  }
1138 }
1139 #endif
mgx::Vector::data
CU_HOST_DEVICE T * data()
Returns a raw pointer on the data.
Definition: Vector.hpp:303
mgx::Matrix::operator<<
friend QTextStream & operator<<(QTextStream &out, const Matrix &mat)
Definition: Matrix.hpp:365
mgx::Matrix::operator*=
CU_HOST_DEVICE Matrix & operator*=(const Matrix &mat)
Definition: Matrix.hpp:322
mgx::Matrix::map
CU_HOST_DEVICE Matrix< nRows, nCols, T > map(T(*fct)(T, T), const Matrix< nRows, nCols, T > &m1, const Matrix< nRows, nCols, T > &m2)
Apply a binary function to each element of the matrix.
Definition: Matrix.hpp:798
mgx::transpose
CU_HOST_DEVICE Matrix< nCols, nRows, T > transpose(const Matrix< nRows, nCols, T > &mat)
Definition: Matrix.hpp:1114
Vector.hpp
mgx::cofactor
CU_HOST_DEVICE T cofactor(const Matrix< nRows, nRows, T > &mat, size_t i, size_t j)
Definition: Matrix.hpp:992
mgx::Matrix::Matrix
CU_HOST_DEVICE Matrix(void)
Create a matrix filled with 0s.
Definition: Matrix.hpp:59
mgx::Matrix::c_data
const CU_HOST_DEVICE T * c_data() const
Returns a constant raw pointer on the data.
Definition: Matrix.hpp:178
mgx::Matrix::operator/
CU_HOST_DEVICE Matrix operator/(const T &scalar) const
Matrix-scalar division.
Definition: Matrix.hpp:252
mgx::Matrix::map
CU_HOST_DEVICE Matrix< nRows, nCols, T > map(T(*fct)(T1), const Matrix< nRows, nCols, T1 > &m)
Apply a unary function to each element of the matrix.
Definition: Matrix.hpp:760
mgx::normsq
CU_HOST_DEVICE T normsq(const Matrix< nRows, nCols, T > &mat)
Definition: Matrix.hpp:1123
mgx::Matrix::Matrix
CU_HOST_DEVICE Matrix(const Vector< nCols, T1 > &r0, const Vector< nCols, T1 > &r1, const Vector< nCols, T1 > &r2)
Initialize a 3 row matrix.
Definition: Matrix.hpp:93
mgx::return
return(nrml ^(c3 *(p2 - p1)+c1 *(p3 - p2)+c2 *(p1 - p3)))/area2
mgx::operator*
const VSMultOp< DistNhbdT, MatrixT, VectorT, ScalarT > operator*(DVector< DistNhbdT, MatrixT, VectorT, ScalarT > &p_vec, ScalarT p_a)
Definition: DistMatrix.hpp:159
StaticAssert.hpp
STATIC_ASSERT
#define STATIC_ASSERT(v)
Definition: StaticAssert.hpp:20
mgx::Matrix::diag
CU_HOST_DEVICE Vector< nRows, T > diag() const
Return the diagonal vector, if the matrix is square.
Definition: Matrix.hpp:560
mgx::Vector::z
CU_HOST_DEVICE void z(const T &v)
Short access to the third element.
Definition: Vector.hpp:739
mgx::Matrix::operator/=
CU_HOST_DEVICE Matrix & operator/=(const T &scalar)
Definition: Matrix.hpp:317
mgx::Matrix::operator*
CU_HOST_DEVICE Vector< nRows, T > operator*(const Vector< nCols, T > &vec) const
Matrix*Column Vector.
Definition: Matrix.hpp:280
mgx::Matrix::operator-
CU_HOST_DEVICE Matrix operator-(void) const
Matrix subtraction.
Definition: Matrix.hpp:196
mgx::Matrix::map
CU_HOST_DEVICE Matrix< nRows, nCols, T > map(T(*fct)(const T &, const T &), const Matrix< nRows, nCols, T > &m1, const Matrix< nRows, nCols, T > &m2)
Apply a binary function to each element of the matrix.
Definition: Matrix.hpp:819
mgx::Matrix::map
CU_HOST_DEVICE Matrix< nRows, nCols, T > map(const T &(*fct)(const T &), const Matrix< nRows, nCols, T > &m)
Apply a unary function to each element of the matrix.
Definition: Matrix.hpp:684
mgx::matrix_minor
CU_HOST_DEVICE T matrix_minor(const Matrix< nRows, nRows, T > &mat, size_t i, size_t j)
Definition: Matrix.hpp:972
mgx::Matrix::Matrix
CU_HOST_DEVICE Matrix(const T &value)
Create a diagonal matrix.
Definition: Matrix.hpp:139
mgx::inverse
CU_HOST_DEVICE Matrix< 1, 1, T > inverse(const Matrix< 1, 1, T > &mat)
Definition: Matrix.hpp:1001
mgx::Matrix::operator+
CU_HOST_DEVICE Matrix operator+(const Matrix &mat) const
Matrix addition.
Definition: Matrix.hpp:210
mgx::Vector::t
CU_HOST_DEVICE void t(const T &v)
Short access to the fourth element.
Definition: Vector.hpp:748
mgx
Distributed matrix library.
Definition: Assert.hpp:26
mgx::Matrix::map
CU_HOST_DEVICE Matrix< nRows, nCols, T > map(const T &(*fct)(const T1 &), const Matrix< nRows, nCols, T1 > &m)
Apply a unary function to each element of the matrix.
Definition: Matrix.hpp:741
mgx::Matrix::rotation
static CU_HOST_DEVICE Matrix< 4, 4, T > rotation(const Vector< 4, T > &direction, T angle)
Creates the 4x4 matrix corresponding to a rotation.
Definition: Matrix.hpp:522
mgx::Matrix::trace
CU_HOST_DEVICE T trace() const
Trace of the matrix.
Definition: Matrix.hpp:541
mgx::Matrix::operator*=
CU_HOST_DEVICE Matrix & operator*=(const T &scalar)
Definition: Matrix.hpp:312
mgx::Matrix::map
CU_HOST_DEVICE Matrix< nRows, nCols, T > map(T(*fct)(const T1 &, const T2 &), const Matrix< nRows, nCols, T1 > &m1, const Matrix< nRows, nCols, T2 > &m2)
Apply a binary function to each element of the matrix.
Definition: Matrix.hpp:882
mgx::Matrix::operator~
CU_HOST_DEVICE Matrix< nCols, nRows, T > operator~()
Transpose the matrix.
Definition: Matrix.hpp:485
mgx::Matrix::Matrix
CU_HOST_DEVICE Matrix(const Vector< nCols, T1 > *vecs)
Fill the matrix with the array of vectors.
Definition: Matrix.hpp:83
mgx::Matrix::operator>>
friend QTextStream & operator>>(QTextStream &in, Matrix &mat)
Definition: Matrix.hpp:376
mgx::Matrix::data
CU_HOST_DEVICE T * data()
Returns a raw pointer on the data.
Definition: Matrix.hpp:188
mgx::Matrix::reference_type
T & reference_type
Definition: Matrix.hpp:45
mgx::Matrix::numcols
static const size_t numcols
Definition: Matrix.hpp:52
mgx::Matrix::operator*
CU_HOST_DEVICE Matrix operator*(const T &scalar) const
Matrix-scalar multiplication.
Definition: Matrix.hpp:238
mgx::Matrix::const_iterator
const typedef T * const_iterator
Definition: Matrix.hpp:50
mgx::Matrix::const_reference_type
const typedef T & const_reference_type
Definition: Matrix.hpp:46
mgx::Matrix::numrows
static const size_t numrows
Definition: Matrix.hpp:53
mgx::Vector::y
CU_HOST_DEVICE void y(const T &v)
Short access to the second element.
Definition: Vector.hpp:730
mgx::det
CU_HOST_DEVICE T det(const Matrix< 1, 1, T > &mat)
Definition: Matrix.hpp:946
mgx::Matrix::Matrix
CU_HOST_DEVICE Matrix(const T *values)
Fill in the matrix with the values.
Definition: Matrix.hpp:127
mgx::Matrix::operator[]
CU_HOST_DEVICE Vector< nCols, T > & operator[](size_t idx)
Returns the nth row.
Definition: Matrix.hpp:410
Util.hpp
mgx::Matrix::operator==
CU_HOST_DEVICE bool operator==(const Matrix &mat) const
Element-wise equality.
Definition: Matrix.hpp:330
mgx::Matrix::size
static CU_HOST_DEVICE Vector< 2, size_t > size()
Returns the size of the matrix.
Definition: Matrix.hpp:152
CudaGlobal.hpp
mgx::Matrix::operator>>
CU_HOST_DEVICE friend std::istream & operator>>(std::istream &in, Matrix &mat)
Definition: Matrix.hpp:397
mgx::Matrix::zero
CU_HOST_DEVICE Matrix & zero(void)
Set the matrix to all zero.
Definition: Matrix.hpp:454
mgx::Matrix::fillArray
CU_HOST_DEVICE void fillArray(T *array)
Definition: Matrix.hpp:551
mgx::Matrix::pointer_type
T * pointer_type
Definition: Matrix.hpp:47
mgx::Matrix::map
CU_HOST_DEVICE Matrix< nRows, nCols, T > map(const T &(*fct)(const T &, const T &), const Matrix< nRows, nCols, T > &m1, const Matrix< nRows, nCols, T > &m2)
Apply a binary function to each element of the matrix.
Definition: Matrix.hpp:840
mgx::Information::out
mgx_EXPORT QTextStream out
mgx::Matrix::rotation
static CU_HOST_DEVICE Matrix< 3, 3, T > rotation(const Vector< 3, T > &direction, T angle)
Creates the 3x3 matrix corresponding to a rotation.
Definition: Matrix.hpp:501
mgx::Vector
Namespace containing all the utility classes.
Definition: Vector.hpp:48
mgx::Matrix::map
CU_HOST_DEVICE Matrix< nRows, nCols, T > map(const T &(*fct)(const T1 &, const T2 &), const Matrix< nRows, nCols, T1 > &m1, const Matrix< nRows, nCols, T2 > &m2)
Apply a binary function to each element of the matrix.
Definition: Matrix.hpp:903
mgx::Matrix::nbRows
static CU_HOST_DEVICE size_t nbRows()
Returns the number of rows of the matrix.
Definition: Matrix.hpp:160
mgx::norm
ScalarT norm(DVector< DistNhbdT, MatrixT, VectorT, ScalarT > &v)
Definition: DistMatrix.hpp:540
mgx::Matrix::operator=
CU_HOST_DEVICE Matrix & operator=(const T &value)
Set the matrix to a diagonal matrix.
Definition: Matrix.hpp:468
mgx::min
CU_HOST_DEVICE T min(const T a, const T b)
Definition: Util.hpp:26
mgx::Matrix::const_pointer_type
const typedef T * const_pointer_type
Definition: Matrix.hpp:48
mgx::Matrix::operator+=
CU_HOST_DEVICE Matrix & operator+=(const Matrix &mat)
Definition: Matrix.hpp:302
CU_HOST_DEVICE
#define CU_HOST_DEVICE
Definition: CudaGlobal.hpp:22
mgx::Vector::c_data
const CU_HOST_DEVICE T * c_data() const
Returns a constant raw pointer on the data.
Definition: Vector.hpp:331
mgx::Matrix::nbColumns
static CU_HOST_DEVICE size_t nbColumns()
Returns the number of columns of the matrix.
Definition: Matrix.hpp:168
mgx::Matrix::operator()
CU_HOST_DEVICE T operator()(size_t i, size_t j) const
Return the value at row i, column j.
Definition: Matrix.hpp:436
mgx::Matrix::iterator
T * iterator
Definition: Matrix.hpp:49
mgx::Matrix::operator!=
CU_HOST_DEVICE bool operator!=(const Matrix &mat) const
Element-wise inequality.
Definition: Matrix.hpp:343
mgx::Matrix::map
CU_HOST_DEVICE Matrix< nRows, nCols, T > map(T(*fct)(const T1 &), const Matrix< nRows, nCols, T1 > &m)
Apply a unary function to each element of the matrix.
Definition: Matrix.hpp:779
mgx::mat2
Matrix< 2, 2, GLfloat > mat2
Definition: Shader.hpp:40
mgx::Matrix::operator=
CU_HOST_DEVICE Matrix & operator=(const Matrix &mat)
Definition: Matrix.hpp:293
mgx::Vector::x
CU_HOST_DEVICE void x(const T &v)
Short access to the first element.
Definition: Vector.hpp:721
mgx::Matrix
Definition: Matrix.hpp:39
mgx::Matrix::operator-=
CU_HOST_DEVICE Matrix & operator-=(const Matrix &mat)
Definition: Matrix.hpp:307
mgx::Matrix::operator-
CU_HOST_DEVICE Matrix operator-(const Matrix &mat) const
Definition: Matrix.hpp:224
mgx::Matrix::operator()
CU_HOST_DEVICE T & operator()(size_t i, size_t j)
Return the value at row i, column j.
Definition: Matrix.hpp:428
mgx::Matrix::operator*
CU_HOST_DEVICE friend Matrix operator*(const T &scalar, const Matrix &mat)
Matrix-scalar multiplication.
Definition: Matrix.hpp:266
mgx::Matrix::map
CU_HOST_DEVICE Matrix< nRows, nCols, T > map(T(*fct)(T), const Matrix< nRows, nCols, T > &m)
Apply a unary function to each element of the matrix.
Definition: Matrix.hpp:703
mgx::Matrix::map
CU_HOST_DEVICE Matrix< nRows, nCols, T > map(T(*fct)(const T &), const Matrix< nRows, nCols, T > &m)
Apply a unary function to each element of the matrix.
Definition: Matrix.hpp:722
mgx::Matrix::identity
static CU_HOST_DEVICE Matrix identity()
Returns an identity matrix.
Definition: Matrix.hpp:444
mgx::Matrix::operator[]
CU_HOST_DEVICE Vector< nCols, T > operator[](size_t idx) const
Returns the nth row.
Definition: Matrix.hpp:420
mgx::Matrix::map
CU_HOST_DEVICE Matrix< nRows, nCols, T > map(T(*fct)(T1, T2), const Matrix< nRows, nCols, T1 > &m1, const Matrix< nRows, nCols, T2 > &m2)
Apply a binary function to each element of the matrix.
Definition: Matrix.hpp:861
mgx::Matrix::operator<
CU_HOST_DEVICE bool operator<(const Matrix &other) const
Comparison operator.
Definition: Matrix.hpp:353
mgx::Matrix::Matrix
CU_HOST_DEVICE Matrix(const T1 *values)
Fill in the matrix with the values.
Definition: Matrix.hpp:111
mgx::Matrix::operator<<
CU_HOST_DEVICE friend std::ostream & operator<<(std::ostream &out, const Matrix &mat)
Definition: Matrix.hpp:385
mgx::Matrix::Matrix
CU_HOST_DEVICE Matrix(const Matrix< nRows, nCols, T1 > &mat)
Copy a matrix.
Definition: Matrix.hpp:71
mgx::Matrix::value_type
T value_type
Definition: Matrix.hpp:44