EMatrix
EMatrix.h
Go to the documentation of this file.
1 /** This file is part of EMatrix, the C++ matrix library distribution.
2  * This project is licensed under the terms of the MIT license. The full text
3  * of the license file can be found in LICENSE.
4  */
5 
6 /// \file
7 
8 #ifndef _EMATRIX_H
9 #define _EMATRIX_H 1
10 
11 #include <cmath>
12 #include <cfloat>
13 #include <cstdlib>
14 #include <cstring>
15 #include <cstdarg>
16 #include <cassert>
17 
18 #include <complex>
19 #include <fstream>
20 #include <iostream>
21 #include <initializer_list>
22 #include <random>
23 #include <limits>
24 
25 
26 #define cout_octave(x) octave((cout),(x),(#x))
27 #define HERE std::cerr << __FILE__ << ':' << __LINE__ << std::endl;
28 
29 namespace ematrix {
30 
31 //! These constants are from the gcc 3.2 <cmath> file (overkill??)
32 namespace cnst {
33 const double PI = 3.1415926535897932384626433832795029L; // pi
34 const double PI_2 = 1.5707963267948966192313216916397514L; // pi/2
35 const double RTD = 180.0L/PI; // Radians to Degrees
36 const double DTR = 1.0L/RTD; // Degrees to Radians
37 
38 /** WGS84 Earth constants:
39  * Titterton, D.H. and Weston, J.L.; <i>Strapdown Inertial Navigation
40  * Technology</i>; Peter Peregrinus Ltd.; London 1997; ISBN 0-86341-260-2; pg. 53.
41  */
42 const double REQ = 6378137.0; // meters
43 const double rEQ = 6356752.3142; // meters
44 const double f = 1.0/298.257223563;
45 const double ECC = 0.0818191908426;
46 const double ECC2 = ECC*ECC;
47 const double WS = 7.292115E-05; // rad/sec
48 
49 const double g = 9.78032534; // m/s^2
50 const double k = 0.00193185; // ??
51 }
52 
53 template < typename tData, size_t tRows, size_t tCols >
54 class Matrix {
55 protected:
56  // Memory allocation method
57  void matalloc (size_t iRowIndex, size_t iColIndex); // Tag:tested
58  // Number of row and columns
59  size_t iRows, iCols;
60  // Storage element, matalloc above assigns an array of pointers to pointers
61  tData *ij[tRows];
62  tData storage[tRows*tCols];
63 public:
64 
65  //! Virtual destructor, no need though
66  virtual ~Matrix (); // Tag:tested
67 
68  /** Default constructor
69  * Usage: Matrix<double,2,3> A;
70  */
71  Matrix (); // Tag:tested
72 
73  /** Copy constructor (not to be confused with the assignment operator)
74  * Usage: Matrix<float,2,3> A;
75  * Matrix<float,2,3> B=A;
76  */
77  inline Matrix (const Matrix< tData, tRows, tCols > & R); // Tag:tested
78 
79  /** Array initialize contructor
80  * Usage: float a[2][3] = {{1.0,2.0,3.0},{4.0,5.0,6.0}};
81  * Matrix<float,2,3> A(&a[0][0]);
82  */
83  Matrix (tData* tArray); // Tag:tested
84 
85  /** STL list initialize contructor (C++11)
86  * Usage: Matrix<double,3,3> A = {1.,2.,3.,0.,0.,0.,0,0,0};
87  */
88  Matrix (const std::initializer_list<tData>& l); // Tag:tested
89 
90  /** Assignment operator (not to be confused with the copy constructor)
91  * Usage: Y = X - Z;
92  */
93  inline const Matrix< tData, tRows, tCols > &operator = (const Matrix< tData, tRows, tCols > & R); // Tag:tested
94 
95  /** STL initializer_list assignment (C++11)
96  * Usage: Matrix<double,3,2> A;
97  * A = {1.1,2.1,3.1,0.0,0.0,0.0};
98  *
99  */
100  inline const Matrix< tData, tRows, tCols > &operator = (const std::initializer_list<tData>& l); // Tag:tested
101 
102  /** Array assignment
103  * Usage: double a[2][3] = {{1.1,2.1,3.1},{4.0,5.0,6.0}};
104  * Matrix<double,2,3> A;
105  * A.load(&a[0][0]);
106  * Warning: Not congnizant of size, can read from unintended memory location;
107  */
108  void load(const tData* tArray); // Tag:tested consider changing to memset
109 
110  /** C like element access (0 to n-1), get and set.
111  * Usage: Matrix<double,3,2> A = {1,2,3,4,5,6};
112  * A[0][0] = 7;
113  * cerr << A[2][1] << endl;
114  * Row operator returns the matrix row corresponding to iRowIndex,
115  * Warning: Does not provide column access safety.
116  */
117  inline tData * operator [] (size_t iRowIndex) const { // Tag:tested
118  assert(iRowIndex<iRows);
119  return ij[iRowIndex];
120  }
121 
122  /** Data access operator for Octave and FORTRAN indexing ... From 1 to n
123  * Note this does not imply FORTRAN memory storage
124  * Usage: Matrix<double,3,2> A = {1,2,3,4,5,6};
125  * A(1,1) = 8;
126  * cerr << A(3,2) << endl;
127  * Note this looks similar to the deprecated memory initialize c_tor that uses va_arg.
128  * but is not the same thing
129  */
130  inline tData & operator () (size_t iRowIndex, size_t iColIndex) const; // Tag:tested
131 
132  /** Vector Data access operator for Octave and FORTRAN indexing ... From 1 to n
133  * Usage: Matrix<double,6,1> V = {1,2,3,4,5,6}; V(1) = 8; cerr << V(6) << endl;
134  * Usage: Matrix<double,1,6> U = {1,2,3,4,5,6}; U(1) = 8; cerr << U(6) << endl;
135  * Note a non 1xn or nx1 matrix will assertion fail. Could not determine a way
136  * to force compile time error.
137  */
138  inline tData & operator () (size_t iIndex) const; // Tag:tested
139 
140  /** Overloaded output stream operator <<
141  * Usage: log_file << A;
142  */
143  template < class tData0, size_t tRows0, size_t tCols0 >
144  friend std::ostream& operator << (std::ostream& s,const Matrix< tData0, tRows0, tCols0 >& A);// Tag:tested
145 
146  /** Octave text output format, complex<double>
147  * Usage: Matrix< complex<double>,3,3 > ZA;
148  * octave(cout,ZA,"ZA") << endl;
149  * Can define: #define cout_octave(x) octave(cout,x,#x) for
150  * cout_octave(ZA) << endl;
151  */
152  // This is broken but works well enough. See:
153  // https://web.mst.edu/~nmjxv3/articles/templates.html
154  template < class tData0, size_t tRows0, size_t tCols0 >
155  friend std::ostream& octave (std::ostream& s, const Matrix< tData0, tRows0, tCols0 >& A, const char* Aname);// { // Tag:tested
156 
157  /** Octave text output format, double
158  * Usage: Matrix<double,2,3> A(&a[0][0]);
159  * octave(cout,A,"A") << endl;
160  * Can define: #define cout_octave(x) octave(cout,x,#x) for
161  * cout_octave(A) << endl;
162  */
163  template < class tData0, size_t tRows0, size_t tCols0 >
164  friend std::ostream& octave (std::ostream& s, const Matrix< std::complex<tData0>, tRows0, tCols0 >& A, const char* Aname);// { // Tag:tested
165 
166  /** Get the storage pointer for the data in the matrix
167  * This is really only here for the friend functions
168  * Try not to use it
169  * Usage: tData* ptr = A.pIJ();
170  */
171  inline tData *pIJ (void) const { // Tag:tested
172  return ij[0];
173  }
174 
175  /** Get the number of rows in a matrix
176  * Usage: size_t i = A.rows();
177  */
178  inline size_t rows (void) const { // Tag:tested
179  return iRows;
180  }
181 
182  /** Get the number of cols in a matrix
183  * Usage: size_t i = A.cols();
184  */
185  inline size_t cols (void) const { // Tag:tested
186  return iCols;
187  }
188 
189  /** Boolean == operator
190  * Usage: if (A == B) ...
191  */
192  bool operator == (Matrix< tData, tRows, tCols > & R); // Tag:tested
193 
194  /** Boolean != operator
195  * Usage: if (A != B) ...
196  */
197  bool operator != (Matrix< tData, tRows, tCols > & R); // Tag:tested
198 
199  /** Unary + operator
200  * Usage: C = (+B); Just returns *this;
201  */
202  Matrix< tData, tRows, tCols > operator + (); // Tag:tested
203 
204  /** Unary - operator
205  * Usage: C = (-A);
206  */
207  Matrix< tData, tRows, tCols > operator - (); // Tag:tested
208 
209  /** Addition operator
210  * Usage: C = A + B;
211  */
212  Matrix< tData, tRows, tCols > operator + (const Matrix< tData, tRows, tCols > & R); // Tag:tested
213 
214  /** Subtaction operator
215  * Usage: C = A - B;
216  */
217  Matrix< tData, tRows, tCols > operator - (const Matrix< tData, tRows, tCols > & R); // Tag:tested
218 
219  /** Scalar multiplication operator
220  * Usage: C = A * scalar;
221  */
222  Matrix< tData, tRows, tCols > operator * (const tData & scalar); // Tag:tested
223 
224  /** Friend scalar multiplication operator
225  * Usage: C = scalar * A;
226  */
227  template< class tData0, size_t tRows0, size_t tCols0 >
228  friend Matrix< tData0, tRows0, tCols0 > operator * (const tData0 & scalar,const Matrix< tData0, tRows0, tCols0 > & R); // Tag:tested
229 
230  /** Scalar division operator
231  * Usage: C = A / scalar;
232  */
233  Matrix< tData, tRows, tCols > operator / (const tData & scalar); // Tag:tested
234 
235  /** Matrix multiplication operator
236  * Usage: C = A * B;
237  */
238  template < size_t tColsR >
240  tData x;
242 
243  for (size_t iIndex=0; iIndex<iRows; iIndex++) {
244  for (size_t jIndex=0; jIndex<R.cols(); jIndex++) {
245  x = tData(0);
246  for (size_t kIndex=0; kIndex<R.rows(); kIndex++) {
247  x += ij[iIndex][kIndex] * R[kIndex][jIndex];
248  }
249  Result[iIndex][jIndex] = x;
250  }
251  }
252  return Result;
253  }
254 
255  /** Array multiplication operator
256  * Usage: C = (A *= B); Must use parentheses
257  * This mimics Octave's A .* B operator and not C's x *= 5 operator
258  */
259  Matrix< tData, tRows, tCols > operator *= (const Matrix< tData, tRows, tCols > & R); // Tag:tested
260 
261  /** Array division operator
262  * Usage: C = (A /= B); Must use parentheses
263  * This mimics Octave's A ./ B operator and not C's x /= 5 operator
264  */
265  Matrix< tData, tRows, tCols > operator /= (const Matrix< tData, tRows, tCols > & R); // Tag:tested
266 
267  /** Concatenate matrices top and botton
268  * Usage: C = (A & B); Must use parenthesis
269  */
270  template < class tData0, size_t tCols0, size_t tRowsT, size_t tRowsB>
272  const Matrix< tData0, tRowsB, tCols0 >& Bottom); // Tag:tested
273 
274  /** Concatenate matrices Left to Right
275  * Usage: C = (A | B); Must use parenthesis
276  */
277  template < class tData0, size_t tRows0, size_t tColsL, size_t tColsR >
279  const Matrix< tData0, tRows0, tColsR >& Right); // Tag:tested
280 
281  /** Set contents to 0x0
282  * Usage: A.zeros();
283  */
284  Matrix< tData, tRows, tCols > zeros( void ); // Tag:tested
285 
286  /** Set contents to tData(1)
287  * Usage: A.ones();
288  */
289  Matrix< tData, tRows, tCols > ones( void ); // Tag:tested
290 
291  /** Set contents to the identity matrix
292  * Usage: A.eye();
293  */
294  Matrix< tData, tRows, tCols > eye( void ); // Tag:tested
295 
296  /** Set all elements of the current matrix random ~N(0,1);
297  * Usage: u.randn();
298  */
299  Matrix< tData, tRows, tCols > randn(void); // Tag:tested
300 
301  // Matrix transpose and complex conjugate transpose
302  // Usage: A_trans = trans(A);
303  template < class tData0, size_t tRows0, size_t tCols0 >
305 
306  template < size_t tRows0, size_t tCols0 >
307  friend Matrix< std::complex<float>, tCols0, tRows0 > trans( const Matrix< std::complex<float>, tRows0, tCols0 >& R ); // Tag:tested
308 
309  template < size_t tRows0, size_t tCols0 >
310  friend Matrix< std::complex<double>, tCols0, tRows0 > trans( const Matrix< std::complex<double>, tRows0, tCols0 >& R ); // Tag:tested
311 
312  /** Matrix diagonal like Octave
313  * This friend function does not modify input contents.
314  * Usage: A_diag = diag(A);
315  */
316  template < class tData0, size_t tRows0 >
317  friend Matrix< tData0, tRows0, 1 > diag( const Matrix< tData0, tRows0, tRows0 >& R ); // Tag:tested
318 
319  template < class tData0, size_t tRows0 >
320  friend Matrix< tData0, tRows0, tRows0 > diag( const Matrix< tData0, tRows0, 1 >& R ); // Tag:tested
321 
322  template < class tData0, size_t tCols0 >
323  friend Matrix< tData0, tCols0, tCols0 > diag( const Matrix< tData0, 1, tCols0 >& R ); // Tag:tested
324 
325 
326  /* Construct a skew symmetric matrix from a 3x1 vector.
327  * w = [wx;wy;wz];
328  * skew(w) = [0.0 -wz +wy
329  * +wz 0.0 -wx
330  * -wy +wx 0.0];
331  * Usage: omega_x = skew(w);
332  */
333  template < class tData0 >
334  friend Matrix< tData0, 3, 3 > skew( const Matrix< tData0, 3, 1 >& R ); // Tag:tested
335 
336  /* Take the cross product of two 3x1 vectors
337  * Usage: a = cross(lsr,Vc);
338  */
339  template < class tData0 >
340  friend Matrix< tData0, 3, 1 > cross( const Matrix< tData0, 3, 1 >& L, const Matrix< tData0, 3, 1 >& R ); // Tag:tested
341 
342  /* Take the dot product of two 3x1 vectors
343  * Usage: (norm(x))^2 = dot(x,x);
344  */
345  template < class tData0, size_t tRows0 >
346  friend tData0 dot( const Matrix< tData0, tRows0, 1 >& L, const Matrix< tData0, tRows0, 1 >& R ); // Tag:tested
347 
348  /** Take the norm of two vectors
349  * Usage: norm_a = a.n();
350  */
351  tData n( void ) const; // Tag:tested
352 
353  /** Take the norm of two vectors
354  * Usage: norm_a = norm(a);
355  */
356  template < class tData0, size_t tRows0 >
357  friend tData0 norm( const Matrix< tData0, tRows0, 1 >& R ); // Tag:tested
358 
359  /** return a unit vector in the direction of V
360  * Usage: u_v = V.u()
361  */
362  Matrix< tData, tRows, 1 > u( void ) const; // Tag:tested
363 
364  /** Matrix inverse, must link with Lapack
365  * Usage: A_inv = inv(A);
366  */
367  template < size_t tRows0 >
369 
370  template < size_t tRows0>
371  friend Matrix< std::complex<float>, tRows0, tRows0 > inv( const Matrix< std::complex<float>, tRows0, tRows0 >& R );
372 
373  template < size_t tRows0 >
375 
376  template < size_t tRows0 >
377  friend Matrix< std::complex<double>, tRows0, tRows0 > inv( const Matrix< std::complex<double>, tRows0, tRows0 >& R );
378 
379  /** Matrix determinent, must link with Lapack
380  * Usage: A_det = det(A);
381  */
382  template < size_t tRows0 >
383  friend float det( const Matrix< float, tRows0, tRows0 >& R );
384 
385  template < size_t tRows0 >
386  friend double det( const Matrix< double, tRows0, tRows0 >& R );
387 
388  template < size_t tRows0 >
389  friend std::complex<float> det( const Matrix< std::complex<float>, tRows0, tRows0 >& R );
390 
391  template < size_t tRows0 >
392  friend std::complex<double> det( const Matrix< std::complex<double>, tRows0, tRows0 >& R );
393 
394  /** Simple Rotation
395  * Usage: Rx = R(angle_rad, 'x');
396  * Note that r_ECEF = trans(R(ws*t, 'z'))*r_ECI if ws*t is a positive rotation
397  * about the z axis. This is backwards on purpose.
398  */
399  template < class tData0 >
400  friend Matrix< tData0, 3, 3 > R(tData0 angle, char axis);
401  //friend Matrix< tData, 3, 3 > TIE(tData t);
402  //friend Matrix< tData, 3, 3 > TEL(tData lat_rad,tData lon_rad);
403  //friend Matrix< tData, 3, 3 > ang_to_tib(tData psi,tData the,tData phi);
404  //friend void uv_to_ae(tData *az,tData *el,Matrix< tData, 3, 1 > &u);
405 
406 }; // EOC: End of Class
407 
408 
409 template < class tData, size_t tRows, size_t tCols >
410 void Matrix< tData, tRows, tCols >::matalloc(size_t iRowIndex, size_t iColIndex) { // Tag:tested
411  ij[0] = &storage[0];
412  for (size_t iIndex = 1; iIndex < iRowIndex; iIndex++)
413  ij[iIndex] = ij[iIndex - 1] + iColIndex;
414 
415  std::memset (storage, 0x0, sizeof(storage));
416 }
417 
418 template < class tData, size_t tRows, size_t tCols >
420  // HERE;
421 }
422 
423 template < class tData, size_t tRows, size_t tCols >
425  : iRows(tRows), iCols(tCols) {
426  matalloc(tRows,tCols);
427 }
428 
429 template < class tData, size_t tRows, size_t tCols >
430 inline Matrix< tData, tRows, tCols >::Matrix(const Matrix & R) // Tag:tested
431  : iRows (R.iRows), iCols (R.iCols) {
432  matalloc(R.iRows, R.iCols);
433  std::memcpy(ij[0], R.ij[0], sizeof (tData) * iRows * iCols);
434 }
435 
436 template < class tData, size_t tRows, size_t tCols >
437 Matrix< tData, tRows, tCols >::Matrix (tData* tArray) // Tag:tested
438  : iRows (tRows), iCols (tCols) {
439  matalloc(tRows,tCols);
440  std::memcpy(storage, tArray, sizeof(storage));
441 }
442 
443 template < class tData, size_t tRows, size_t tCols >
444 Matrix< tData, tRows, tCols >::Matrix (const std::initializer_list<tData>& l) // Tag:tested
445  : iRows(tRows), iCols(tCols) {
446  assert( iRows*iCols == l.size() );
447  matalloc( iRows, iCols);
448  for(size_t i = 0; i<(iRows*iCols); i++) {
449  ij[0][i] = *(l.begin() + i);
450  }
451 }
452 
453 template < class tData, size_t tRows, size_t tCols >
455  assert((iRows == R.iRows) && (iCols == R.iCols));
456  if( this != &R ) {
457  std::memcpy (ij[0], R.ij[0], sizeof(tData)*iRows*iCols);
458  }
459  return *this;
460 }
461 
462 template < class tData, size_t tRows, size_t tCols >
463 inline const Matrix< tData, tRows, tCols>& Matrix< tData,tRows,tCols>::operator = (const std::initializer_list<tData>& l) { // Tag:tested
464  assert( iRows*iCols == l.size() );
465  for(size_t i = 0; i<(iRows*iCols); i++) {
466  ij[0][i] = *(l.begin() + i);
467  }
468  return *this;
469 }
470 
471 template < class tData, size_t tRows, size_t tCols >
472 void Matrix< tData, tRows, tCols >::load(const tData* tArray) { // Tag:tested
473  std::memcpy(storage, tArray, sizeof(storage));
474 }
475 
476 template < class tData, size_t tRows, size_t tCols >
477 inline tData & Matrix< tData, tRows, tCols >::operator () (size_t iRowIndex, size_t iColIndex) const { // Tag:tested
478  iRowIndex--;
479  iColIndex--;
480  assert(0<=iRowIndex && iRowIndex<iRows);
481  assert(0<=iColIndex && iColIndex<iCols);
482  return ij[iRowIndex][iColIndex];
483 }
484 
485 template < class tData, size_t tRows, size_t tCols >
486 inline tData & Matrix< tData, tRows, tCols >::operator () (size_t iIndex) const { // Tag:tested
487  iIndex--;
488  assert(iRows==1 || iCols==1);
489  if(iCols == 1) {
490  assert(0<=iIndex && iIndex<iRows);
491  return ij[iIndex][0];
492  } else { // (iRows == 1)
493  assert(0<=iIndex && iIndex<iCols);
494  return ij[0][iIndex];
495  }
496 }
497 
498 template < class tData, size_t tRows, size_t tCols >
499 std::ostream& operator << (std::ostream& s,const Matrix< tData, tRows, tCols >& A) { // Tag:tested
500  // Sets new precision, returns old. Should figure out how to modify
501  // without code change here. Switch to exponential as well.
502  std::streamsize old_precision = s.precision(8);
503 
504  //s.setf( std::ios::fixed, std:: ios::floatfield ); // floatfield set to fixed
505 
506  //s << "Address: 0x" << (&A) << std::endl;
507  for (size_t i=0; i<A.iRows; i++) {
508  for (size_t j=0; j<A.iCols; j++) {
509  //s.width(25);
510  s << (A[i][j]) << '\t';
511  }
512  s << std::endl;
513  }
514  s.flush();
515  s.precision(old_precision);
516  return s;
517 }
518 
519 template < class tData0, size_t tRows0, size_t tCols0 >
520 std::ostream& octave (std::ostream& s, const Matrix< std::complex< tData0 >, tRows0, tCols0 >& A, const char* Aname) { // Tag:tested
521  s << "# name: " << Aname << std::endl;
522  s << "# type: complex matrix" << std::endl;
523  s << "# rows: " << A.iRows << std::endl;
524  s << "# columns: " << A.iCols << std::endl;
525  s << A;
526  s.flush();
527  return s;
528 }
529 
530 template < class tData0, size_t tRows0, size_t tCols0 >
531 std::ostream& octave (std::ostream& s,const Matrix< tData0, tRows0, tCols0 >& A, const char* Aname) { // Tag:tested
532  s << "# name: " << Aname << std::endl;
533  s << "# type: matrix" << std::endl;
534  s << "# rows: " << A.iRows << std::endl;
535  s << "# columns: " << A.iCols << std::endl;
536  s << A;
537  s.flush();
538  return s;
539 }
540 
541 template < class tData, size_t tRows, size_t tCols >
543  tData * pLeft = ij[0];
544  tData * pRight = R.ij[0];
545 
546  for (size_t iIndex = 0; iIndex < iRows * iCols; iIndex++)
547  if((*pLeft++) != tData(*pRight++)) {
548  return false;
549  }
550 
551  return true;
552 }
553 
554 template < class tData, size_t tRows, size_t tCols >
556  return !(*this == R);
557 }
558 
559 template < class tData, size_t tRows, size_t tCols >
561  return *this;
562 }
563 
564 template < class tData, size_t tRows, size_t tCols >
567  tData * pLeft = ij[0];
568  tData * pResult = Result.ij[0];
569  for (size_t iIndex = 0; iIndex < iRows * iCols; iIndex++)
570  (*pResult++) = (-(*pLeft++));
571  return Result;
572 }
573 
574 template < class tData, size_t tRows, size_t tCols >
577  tData * pLeft = ij[0];
578  tData * pRight = R.ij[0];
579  tData * pResult = Result.ij[0];
580 
581  for (size_t iIndex = 0; iIndex < iRows * iCols; iIndex++)
582  (*pResult++) = (*pLeft++) + (*pRight++);
583 
584  return Result;
585 }
586 
587 template < class tData, size_t tRows, size_t tCols >
590  tData * pLeft = ij[0];
591  tData * pRight = R.ij[0];
592  tData * pResult = Result.ij[0];
593 
594  for (size_t iIndex = 0; iIndex < iRows * iCols; iIndex++)
595  (*pResult++) = (*pLeft++) - (*pRight++);
596  return Result;
597 }
598 
599 template < class tData, size_t tRows, size_t tCols >
601  Matrix< tData, tRows, tCols > Result = (*this);
602  tData * pResult = Result.ij[0];
603  for (size_t iIndex = 0; iIndex < iRows * iCols; iIndex++)
604  (*pResult++) *= scalar;
605  return Result;
606 }
607 
608 template < class tData, size_t tRows, size_t tCols > // friend
609 Matrix< tData, tRows, tCols > operator * (const tData & scalar,const Matrix< tData, tRows, tCols > & R) { // Tag:tested
610  size_t iRows = R.iRows;
611  size_t iCols = R.iCols;
613  tData * pResult = Result.ij[0];
614  for (size_t iIndex = 0; iIndex < iRows * iCols; iIndex++)
615  (*pResult++) *= scalar;
616  return Result;
617 }
618 
619 template < class tData, size_t tRows, size_t tCols >
621  assert(scalar);
622  Matrix< tData, tRows, tCols > Result = (*this);
623  tData * pResult = Result.ij[0];
624  for (size_t iIndex = 0; iIndex < iRows * iCols; iIndex++)
625  (*pResult++) /= scalar;
626  return Result;
627 }
628 
629 template < class tData, size_t tRows, size_t tCols >
632  tData * pLeft = ij[0];
633  tData * pRight = R.ij[0];
634  tData * pResult = Result.ij[0];
635  for (size_t iIndex = 0; iIndex < iRows*iCols; iIndex++)
636  (*pResult++) = (*pLeft++) * tData(*pRight++);
637  return Result;
638 }
639 
640 template < class tData, size_t tRows, size_t tCols >
643  tData * pLeft = ij[0];
644  tData * pRight = R.ij[0];
645  tData * pResult = Result.ij[0];
646  for (size_t iIndex = 0; iIndex < iRows * iCols; iIndex++) {
647  assert(*pRight != 0.0);
648  (*pResult++) = (*pLeft++) / tData(*pRight++);
649  }
650  return Result;
651 }
652 
653 template < class tData0, size_t tCols0, size_t tRowsT, size_t tRowsB>
655  const Matrix< tData0, tRowsB, tCols0 >& Bottom) { // Tag:tested
657 
658  size_t i,ii,j;
659  for(i=0; i<Top.iRows; i++)
660  for(j=0; j<Top.iCols; j++)
661  Result[i][j] = Top[i][j];
662 
663  for(i=Top.iRows,ii=0; i<(Top.iRows+Bottom.iRows); i++,ii++) {
664  for(j=0; j<Top.iCols; j++) {
665  Result[i][j] = Bottom[ii][j];
666  }
667  }
668  return Result;
669 }
670 
671 template < class tData0, size_t tRows0, size_t tColsL, size_t tColsR >
673  const Matrix< tData0, tRows0, tColsR >& Right) { // Tag:tested
675  size_t i,j,jj;
676  for(i=0; i<Left.iRows; i++)
677  for(j=0; j<Left.iCols; j++)
678  Result[i][j] = Left[i][j];
679  for(i=0; i<Left.iRows; i++) {
680  for(j=Left.iCols,jj=0; j<(Left.iCols+Right.iCols); j++,jj++) {
681  Result[i][j] = Right[i][jj];
682  }
683  }
684  return Result;
685 }
686 
687 template < class tData, size_t tRows, size_t tCols >
689  std::memset (ij[0], 0x0, sizeof(tData) * iRows * iCols);
690  return *this;
691 }
692 
693 template < class tData, size_t tRows, size_t tCols >
695  tData * pThis = ij[0];
696  for (size_t iIndex = 0; iIndex < iRows * iCols; iIndex++)
697  (*pThis++) = tData(1);
698  return *this;
699 }
700 
701 template < class tData, size_t tRows, size_t tCols >
703  assert(tRows==tCols);
704  std::memset (ij[0], 0x0, sizeof(tData) * iRows * iCols);
705 
706  tData * pThis = ij[0];
707  for (size_t iIndex = 0; iIndex < iRows; iIndex++, pThis+=iCols)
708  (*pThis++) = tData(1);
709  return *this;
710 }
711 
712 template < class tData, size_t tRows, size_t tCols >
714  std::default_random_engine generator;
715  std::normal_distribution<tData> distribution;
716  for (size_t iIndex = 0; iIndex < iRows*iCols; iIndex++) {
717  ij[0][iIndex] = distribution(generator);
718  }
719  return *this;
720 }
721 
722 template < class tData, size_t tRows, size_t tCols >
725 
726  for (size_t iIndex = 0; iIndex < tCols; iIndex++)
727  for (size_t jIndex = 0; jIndex < tRows; jIndex++)
728  Result[iIndex][jIndex] = R[jIndex][iIndex];
729 
730  return Result;
731 }
732 
733 template < size_t tRows, size_t tCols >
734 Matrix< std::complex<float>, tCols, tRows > trans( const Matrix< std::complex<float>, tRows, tCols >& R ) { // Tag:tested
735  Matrix< std::complex<float>, tCols, tRows > Result;
736 
737  for (size_t iIndex = 0; iIndex < tCols; iIndex++)
738  for (size_t jIndex = 0; jIndex < tRows; jIndex++)
739  Result[iIndex][jIndex] = std::conj(R[jIndex][iIndex]);
740 
741  return Result;
742 }
743 
744 template < size_t tRows, size_t tCols >
745 Matrix< std::complex<double>, tCols, tRows > trans( const Matrix< std::complex<double>, tRows, tCols >& R ) { // Tag:tested
746  Matrix< std::complex<double>, tCols, tRows > Result;
747 
748  for (size_t iIndex = 0; iIndex < tCols; iIndex++)
749  for (size_t jIndex = 0; jIndex < tRows; jIndex++)
750  Result[iIndex][jIndex] = std::conj(R[jIndex][iIndex]);
751 
752  return Result;
753 }
754 
755 template < class tData, size_t tRows >
758  for (size_t iIndex = 0; iIndex < tRows; iIndex++ ) {
759  Result[iIndex][0] = R[iIndex][iIndex];
760  }
761  return Result;
762 }
763 
764 template < class tData, size_t tRows >
767  for (size_t iIndex = 0; iIndex < tRows; iIndex++ ) {
768  Result[iIndex][iIndex] = R[iIndex][0];
769  }
770  return Result;
771 }
772 
773 template < class tData, size_t tCols >
776  for (size_t iIndex = 0; iIndex < tCols; iIndex++ ) {
777  Result[iIndex][iIndex] = R[0][iIndex];
778  }
779  return Result;
780 }
781 
782 template < class tData >
784  tData * pR = R.pIJ();
785  tData z = tData(0);
786  Matrix< tData, 3, 3 > Result = {
787  z,-pR[2], pR[1],
788  pR[2], z,-pR[0],
789  -pR[1], pR[0], z
790  };
791 
792  return Result;
793 }
794 
795 template < class tData >
797  return skew(L)*R;
798 }
799 
800 template < class tData, size_t tRows >
801 tData dot( const Matrix< tData, tRows, 1 >& L, const Matrix< tData, tRows, 1 >& R ) { // Tag:tested
802  tData Result = tData(0);
803 
804  tData * pR = R.pIJ();
805  tData * pL = L.pIJ();
806 
807  for (size_t iIndex = 0; iIndex < tRows; iIndex++) {
808  Result += (*pR++)*(*pL++);
809  }
810  return Result;
811 }
812 
813 template < class tData, size_t tRows, size_t tCols >
814 tData Matrix< tData, tRows, tCols >::n( void ) const { // Tag:tested
815  return std::sqrt(dot(*this,*this));
816 }
817 
818 template < class tData, size_t tRows >
819 tData norm( const Matrix< tData, tRows, 1 >& R ) { // Tag:tested
820  return std::sqrt(dot(R,R));
821 }
822 
823 template < class tData, size_t tRows, size_t tCols >
825  tData den = n();
826  assert(den != 0.0);
827  Matrix< tData, tRows, 1 > Result = *this;
828  return Result*(1.0/den);
829 }
830 
831 extern "C" void sgesv_( const int &n, const int &nrhs, float *A,
832  const int &lda, int* ipiv, float *B, const int &ldb, int *info);
833 extern "C" void cgesv_( const int &n, const int &nrhs, std::complex<float> *A,
834  const int &lda, int* ipiv, std::complex<float> *B, const int &ldb, int *info);
835 extern "C" void dgesv_( const int &n, const int &nrhs, double *A,
836  const int &lda, int* ipiv, double *B, const int &ldb, int *info);
837 extern "C" void zgesv_( const int &n, const int &nrhs, std::complex<double> *A,
838  const int &lda, int* ipiv, std::complex<double> *B, const int &ldb, int *info);
839 
840 template < size_t tRows >
842  int n = tRows;
844  int ipiv[tRows] = {0};
846  Result.eye();
847  int info = 0;
848 
849  sgesv_(n, n, a.pIJ(), n, ipiv, Result.pIJ(), n, &info);
850 
851  if(info != 0) {
852  std::cerr << "?gesv returned error: " << info << std::endl;
853  abort();
854  }
855 
856  return Result;
857 }
858 
859 template < size_t tRows >
860 Matrix< std::complex<float>, tRows, tRows > inv( const Matrix< std::complex<float>, tRows, tRows >& R ) {
861  int n = tRows;
862  Matrix< std::complex<float>, tRows, tRows > a = R;
863  int ipiv[tRows] = {0};
864  Matrix< std::complex<float>, tRows, tRows > Result;
865  Result.eye();
866  int info = 0;
867 
868  cgesv_(n, n, a.pIJ(), n, ipiv, Result.pIJ(), n, &info);
869 
870  if(info != 0) {
871  std::cerr << "?gesv returned error: " << info << std::endl;
872  abort();
873  }
874 
875  return Result;
876 }
877 
878 template < size_t tRows >
880  int n = tRows;
882  int ipiv[tRows] = {0};
884  Result.eye();
885  int info = 0;
886 
887  dgesv_(n, n, a.pIJ(), n, ipiv, Result.pIJ(), n, &info);
888 
889  if(info != 0) {
890  std::cerr << "?gesv returned error: " << info << std::endl;
891  abort();
892  }
893 
894  return Result;
895 }
896 
897 template < size_t tRows >
898 Matrix< std::complex<double>, tRows, tRows > inv( const Matrix< std::complex<double>, tRows, tRows >& R ) {
899  int n = tRows;
900  Matrix< std::complex<double>, tRows, tRows > a = R;
901  int ipiv[tRows] = {0};
902  Matrix< std::complex<double>, tRows, tRows > Result;
903  Result.eye();
904  int info = 0;
905 
906  zgesv_(n, n, a.pIJ(), n, ipiv, Result.pIJ(), n, &info);
907 
908  if(info != 0) {
909  std::cerr << "?gesv returned error: " << info << std::endl;
910  abort();
911  }
912 
913  return Result;
914 }
915 
916 extern "C" void sgetrf_(const int &m, const int &n, float *A,
917  const int &lda, int *ipiv, int *info);
918 extern "C" void cgetrf_(const int &m, const int &n, std::complex<float> *A,
919  const int &lda, int *ipiv, int *info);
920 extern "C" void dgetrf_(const int &m, const int &n, double *A,
921  const int &lda, int *ipiv, int *info);
922 extern "C" void zgetrf_(const int &m, const int &n, std::complex<double> *A,
923  const int &lda, int *ipiv, int *info);
924 
925 template < size_t tRows >
927  int n = tRows;
929  int ipiv[tRows] = {0};
930  int info = 0;
931  float result = 1.0;
932 
933  sgetrf_(n,n,a.pIJ(),n,&ipiv[0],&info);
934 
935  if (info == 0) {
936  for(int i=0; i<n; i++) {
937  if(ipiv[i] != (i+1)) result *= -a[i][i]; // i+1 for fortran
938  else result *= +a[i][i];
939  }
940  } else {
941  std::cerr << "?getrf returned error: " << info << std::endl;
942  abort();
943  }
944  return result;
945 }
946 
947 template < size_t tRows >
948 std::complex<float> det( const Matrix< std::complex<float>, tRows, tRows >& R ) {
949  int n = tRows;
950  Matrix< std::complex<float>, tRows, tRows > a = R;
951  int ipiv[tRows] = {0};
952  int info = 0;
953  std::complex<float> result(1.0,0.0);
954 
955  cgetrf_(n,n,a.pIJ(),n,&ipiv[0],&info);
956 
957  if (info == 0) {
958  for(int i=0; i<n; i++) {
959  if(ipiv[i] != (i+1)) result *= -a[i][i]; // i+1 for fortran
960  else result *= +a[i][i];
961  }
962  } else {
963  std::cerr << "?getrf returned error: " << info << std::endl;
964  abort();
965  }
966  return result;
967 }
968 
969 template < size_t tRows >
971  int n = tRows;
973  int ipiv[tRows] = {0};
974  int info = 0;
975  double result = 1.0;
976 
977  dgetrf_(n,n,a.pIJ(),n,&ipiv[0],&info);
978 
979  if (info == 0) {
980  for(int i=0; i<n; i++) {
981  if(ipiv[i] != (i+1)) result *= -a[i][i]; // i+1 for fortran
982  else result *= +a[i][i];
983  }
984  } else {
985  std::cerr << "?getrf returned error: " << info << std::endl;
986  abort();
987  }
988  return result;
989 }
990 
991 template < size_t tRows >
992 std::complex<double> det( const Matrix< std::complex<double>, tRows, tRows >& R ) {
993  int n = tRows;
994  Matrix< std::complex<double>, tRows, tRows > a = R;
995  int ipiv[tRows] = {0};
996  int info = 0;
997  std::complex<double> result(1.0,0.0);
998 
999  zgetrf_(n,n,a.pIJ(),n,&ipiv[0],&info);
1000 
1001  if (info == 0) {
1002  for(int i=0; i<n; i++) {
1003  if(ipiv[i] != (i+1)) result *= -a[i][i]; // i+1 for fortran
1004  else result *= +a[i][i];
1005  }
1006  } else {
1007  std::cerr << "?getrf returned error: " << info << std::endl;
1008  abort();
1009  }
1010  return result;
1011 }
1012 
1013 template < class tData >
1014 Matrix< tData, 3, 3 > R(tData angle, char axis) {
1015  assert(axis == 'x' || axis == 'y' || axis == 'z');
1016  Matrix< tData, 3, 3 > result;
1017 
1018  if(axis == 'x') {
1019  result = {1.0,0.0,0.0,0.0,cos(angle),-sin(angle),0.0,sin(angle), cos(angle)};
1020  } else if(axis=='y') {
1021  result = {cos(angle),0.0,sin(angle),0.0,1.0,0.0,-sin(angle),0.0,cos(angle)};
1022  } else if(axis=='z') {
1023  result = {cos(angle),-sin(angle),0.0,sin(angle),cos(angle),0.0,0.0,0.0,1.0};
1024  } else {
1025  abort();
1026  }
1027  return result;
1028 }
1029 
1030 template < class tData >
1032  return R(cnst::WS*t,'z');
1033 }
1034 
1035 template < class tData >
1036 Matrix< tData, 3, 3 > TEL(tData lat_rad,tData lon_rad) {
1037  return (R(lon_rad,'z')*R(-lat_rad-cnst::PI_2,'y'));
1038 }
1039 
1040 template < class tData >
1041 Matrix< tData, 3, 3 > ang_to_tib(tData psi,tData the,tData phi) {
1042  tData a_r_psi[] = {cos(psi),-sin(psi), 0,
1043  sin(psi), cos(psi), 0,
1044  0, 0, 1
1045  };
1046  Matrix< tData, 3, 3 > r_psi(a_r_psi);
1047 
1048  tData a_r_the[] = {cos(the), 0, sin(the),
1049  0, 1, 0,
1050  -sin(the), 0, cos(the)
1051  };
1052  Matrix< tData, 3, 3 > r_the(a_r_the);
1053 
1054  tData a_r_phi[] = { 1, 0, 0,
1055  0, cos(phi),-sin(phi),
1056  0, sin(phi), cos(phi)
1057  };
1058  Matrix< tData, 3, 3 > r_phi(a_r_phi);
1059 
1060  return r_psi*r_the*r_phi;
1061 }
1062 
1063 template < class tData >
1064 void uv_to_ae(tData *az,tData *el,Matrix< tData, 3, 1 > &u) {
1065  if(norm(u) != 1.0) u=u*(1.0/norm(u));
1066  *az=atan2(u(2),u(1));
1067  *el=asin(-u(3));
1068  return;
1069 }
1070 
1071 } // from namespace
1072 
1073 #endif // from _EMATRIX_H
void uv_to_ae(tData *az, tData *el, Matrix< tData, 3, 1 > &u)
Definition: EMatrix.h:1064
size_t iCols
Definition: EMatrix.h:59
Matrix< tData, tRows, tCols > ones(void)
Definition: EMatrix.h:694
friend Matrix< tData0, 3, 3 > R(tData0 angle, char axis)
size_t iRows
Definition: EMatrix.h:59
std::ostream & octave(std::ostream &s, const Matrix< std::complex< tData0 >, tRows0, tCols0 > &A, const char *Aname)
Definition: EMatrix.h:520
bool operator==(Matrix< tData, tRows, tCols > &R)
Definition: EMatrix.h:542
Matrix< tData, tRows, tCols > operator*(const tData &scalar, const Matrix< tData, tRows, tCols > &R)
Definition: EMatrix.h:609
tData dot(const Matrix< tData, tRows, 1 > &L, const Matrix< tData, tRows, 1 > &R)
Definition: EMatrix.h:801
const double g
Definition: EMatrix.h:49
const double WS
Definition: EMatrix.h:47
Matrix< tData, tRows, tCols > operator*=(const Matrix< tData, tRows, tCols > &R)
Definition: EMatrix.h:630
bool operator!=(Matrix< tData, tRows, tCols > &R)
Definition: EMatrix.h:555
tData n(void) const
Definition: EMatrix.h:814
void zgesv_(const int &n, const int &nrhs, std::complex< double > *A, const int &lda, int *ipiv, std::complex< double > *B, const int &ldb, int *info)
void cgesv_(const int &n, const int &nrhs, std::complex< float > *A, const int &lda, int *ipiv, std::complex< float > *B, const int &ldb, int *info)
Matrix< tData, 3, 3 > TIE(tData t)
Definition: EMatrix.h:1031
Matrix< tData, 3, 3 > R(tData angle, char axis)
Definition: EMatrix.h:1014
Matrix< tData, tRows, tCols > operator+()
Definition: EMatrix.h:560
void dgetrf_(const int &m, const int &n, double *A, const int &lda, int *ipiv, int *info)
Matrix< tData, tRows, tCols > operator/(const tData &scalar)
Definition: EMatrix.h:620
Matrix< tData, tRows, tCols > eye(void)
Definition: EMatrix.h:702
const double f
Definition: EMatrix.h:44
size_t cols(void) const
Definition: EMatrix.h:185
tData * ij[tRows]
Definition: EMatrix.h:61
const double rEQ
Definition: EMatrix.h:43
Matrix< tData, 3, 1 > cross(const Matrix< tData, 3, 1 > &L, const Matrix< tData, 3, 1 > &R)
Definition: EMatrix.h:796
Matrix< tData, 3, 3 > TEL(tData lat_rad, tData lon_rad)
Definition: EMatrix.h:1036
Matrix< tData, tRows, 1 > u(void) const
Definition: EMatrix.h:824
friend tData0 dot(const Matrix< tData0, tRows0, 1 > &L, const Matrix< tData0, tRows0, 1 > &R)
float det(const Matrix< float, tRows, tRows > &R)
Definition: EMatrix.h:926
Matrix< tData, tRows, tCols > zeros(void)
Definition: EMatrix.h:688
friend Matrix< tData0, 3, 1 > cross(const Matrix< tData0, 3, 1 > &L, const Matrix< tData0, 3, 1 > &R)
friend Matrix< tData0, tCols0, tRows0 > trans(const Matrix< tData0, tRows0, tCols0 > &R)
const double PI_2
Definition: EMatrix.h:34
Matrix< float, tRows, tRows > inv(const Matrix< float, tRows, tRows > &R)
Definition: EMatrix.h:841
friend Matrix< float, tRows0, tRows0 > inv(const Matrix< float, tRows0, tRows0 > &R)
void zgetrf_(const int &m, const int &n, std::complex< double > *A, const int &lda, int *ipiv, int *info)
Matrix< tData, tRows, 1 > diag(const Matrix< tData, tRows, tRows > &R)
Definition: EMatrix.h:756
Matrix< tData, 3, 3 > skew(const Matrix< tData, 3, 1 > &R)
Definition: EMatrix.h:783
Matrix< tData, 3, 3 > ang_to_tib(tData psi, tData the, tData phi)
Definition: EMatrix.h:1041
tData & operator()(size_t iRowIndex, size_t iColIndex) const
Definition: EMatrix.h:477
Matrix< tData, tRows, tCols > randn(void)
Definition: EMatrix.h:713
const double ECC
Definition: EMatrix.h:45
Matrix< tData0, tRowsT+tRowsB, tCols0 > operator&(const Matrix< tData0, tRowsT, tCols0 > &Top, const Matrix< tData0, tRowsB, tCols0 > &Bottom)
Definition: EMatrix.h:654
tData * pIJ(void) const
Definition: EMatrix.h:171
size_t rows(void) const
Definition: EMatrix.h:178
friend Matrix< tData0, 3, 3 > skew(const Matrix< tData0, 3, 1 > &R)
void sgesv_(const int &n, const int &nrhs, float *A, const int &lda, int *ipiv, float *B, const int &ldb, int *info)
const Matrix< tData, tRows, tCols > & operator=(const Matrix< tData, tRows, tCols > &R)
Definition: EMatrix.h:454
tData storage[tRows *tCols]
Definition: EMatrix.h:62
friend std::ostream & octave(std::ostream &s, const Matrix< tData0, tRows0, tCols0 > &A, const char *Aname)
Definition: EMatrix.h:531
friend Matrix< tData0, tRows0, 1 > diag(const Matrix< tData0, tRows0, tRows0 > &R)
const double ECC2
Definition: EMatrix.h:46
Matrix< tData, tRows, tCols > operator*(const tData &scalar)
Definition: EMatrix.h:600
void dgesv_(const int &n, const int &nrhs, double *A, const int &lda, int *ipiv, double *B, const int &ldb, int *info)
const double REQ
Definition: EMatrix.h:42
void load(const tData *tArray)
Definition: EMatrix.h:472
Matrix< tData, tRows, tCols > operator-()
Definition: EMatrix.h:565
Matrix< tData0, tRows0, tColsL+tColsR > operator|(const Matrix< tData0, tRows0, tColsL > &Left, const Matrix< tData0, tRows0, tColsR > &Right)
Definition: EMatrix.h:672
const double k
Definition: EMatrix.h:50
Matrix< tData, tRows, tCols > operator/=(const Matrix< tData, tRows, tCols > &R)
Definition: EMatrix.h:641
friend Matrix< tData0, tRows0, tColsL+tColsR > operator|(const Matrix< tData0, tRows0, tColsL > &Left, const Matrix< tData0, tRows0, tColsR > &Right)
Definition: EMatrix.h:672
void cgetrf_(const int &m, const int &n, std::complex< float > *A, const int &lda, int *ipiv, int *info)
const double DTR
Definition: EMatrix.h:36
void matalloc(size_t iRowIndex, size_t iColIndex)
Definition: EMatrix.h:410
friend tData0 norm(const Matrix< tData0, tRows0, 1 > &R)
virtual ~Matrix()
Virtual destructor, no need though.
Definition: EMatrix.h:419
Matrix< tData, tCols, tRows > trans(const Matrix< tData, tRows, tCols > &R)
Definition: EMatrix.h:723
const double RTD
Definition: EMatrix.h:35
friend float det(const Matrix< float, tRows0, tRows0 > &R)
void sgetrf_(const int &m, const int &n, float *A, const int &lda, int *ipiv, int *info)
tData norm(const Matrix< tData, tRows, 1 > &R)
Definition: EMatrix.h:819
friend Matrix< tData0, tRowsT+tRowsB, tCols0 > operator&(const Matrix< tData0, tRowsT, tCols0 > &Top, const Matrix< tData0, tRowsB, tCols0 > &Bottom)
Definition: EMatrix.h:654
const double PI
Definition: EMatrix.h:33