CEBL  2.1
cppR_math.hpp
Go to the documentation of this file.
1 /*
2 * CEBL : CSU EEG Brain-Computer Interface Lab
3 *
4 * Author: Jeshua Bratman - jeshuabratman@gmail.com
5 *
6 * This file is part of CEBL.
7 *
8 * CEBL is free software; you can redistribute it and/or modify it.
9 * We only ask that if you use our code that you cite the source in
10 * your project or publication.
11 *
12 * EEG Group (www.cs.colostate.edu/eeg)
13 * Department of Computer Science
14 * Colorado State University
15 *
16 */
17 
18 
19 //--------------------------------------------------
20 // cppR_math.hpp
21 // Matrix and Vector utilities
22 
23 #ifndef CPPR_MATH_H //make sure this only gets included once
24 #define CPPR_MATH_H
25 
26 #include <cppR/cppR_includes.hpp>
28 
29 using std::complex;
30 
31 //--------------------------------------------------
32 namespace cppR
33 {
34  //--------------------------------------------------
35  //cppR math DATA STRUCTURES
36 
40  template < typename T >
41  struct EigStruct
42  {
43  ublas::matrix<T> vectors;
44  ublas::vector<T> values;
45  };
46 
49  template < typename T >
50  struct SvdStruct
51  {
52  ublas::vector<T> d;
53  ublas::matrix<T> u;
54  ublas::matrix<T> v;
55  };
56 
57  //==============================
65  template < typename T >
66  ublas::matrix<T> compProd(const ublas::matrix<T> & m,
67  const ublas::matrix<T> & n)
68  {
69  cppR_assert(m.size1() == n.size1() && m.size2() == n.size2(),
70  "in compProd: matrices must be of the same size.");
71 
72  ublas::matrix<T> ret(m.size1(),m.size2());
73 
74  for(unsigned int row=0; row<m.size1(); row++)
75  for(unsigned int col=0; col<m.size2(); col++)
76  ret(row,col) = m(row,col) * n(row,col);
77 
78  return ret;
79  }
80 
81  //==============================
82  //compDiv
83 
84  //==============================
92  template < typename T >
93  ublas::matrix<T>
94  compDiv(const ublas::matrix<T> & m, const ublas::matrix<T> & n)
95  {
96  assert(m.size1() == n.size1() && m.size2() == n.size2());
97 
98  ublas::matrix<T> ret(m.size1(),m.size2());
99 
100  for(unsigned row=0; row<m.size1(); row++)
101  for(unsigned col=0; col<m.size2(); col++)
102  ret(row,col) = m(row,col) / n(row,col);
103 
104  return ret;
105  }
106 
107 
108  //==============================
118  template < typename T >
119  ublas::matrix<T> solve(const ublas::matrix<T> & m1)
120  {
121  int size = m1.size1();
122  ublas::matrix<T, ublas::column_major> M(size,size);
123  ublas::matrix<T, ublas::column_major> temp(size,size);
124  M = m1;
125  temp = diag(1,size);
126  lapack::gesv(M,temp);
127 
128  //return a row_major matrix
129  ublas::matrix<T> ret = temp;
130  return ret;
131 
132  };
133 
134  /* Computes solution to set of real linear equations in matrix form.
135  Uses m2 as right-hand side of equations.
136 
137  cppR::solve uses lapack routine gesv to perform the inverse.
138  \param m1
139  \param m2
140 
141  \return
142  */
143  template < typename T >
144  ublas::matrix<T> solve(const ublas::matrix<T> & m1,
145  const ublas::matrix<T> & m2)
146  {
147  int size = m1.size1();
148  ublas::matrix<T, ublas::column_major> M(size,size);
149  ublas::matrix<T, ublas::column_major> temp(size,size);
150  M = m1;
151  temp = m2;
152  lapack::gesv(M,temp);
153 
154  //return a row_major matrix
155  ublas::matrix<T> ret = temp;
156  return ret;
157 
158  };
159 
160  //==============================
167  template < typename T >
168  EigStruct<T> eigen(const ublas::matrix<T> & m1)
169  {
170  int size = m1.size1();
171  ublas::matrix<T, ublas::column_major> M(size,size);
172  ublas::matrix<double, ublas::column_major> B(size,size);
173  ublas::vector<double> D(size);
174  M = m1;
175  lapack::syev('V','U',M,D,lapack::optimal_workspace());
176 
177  EigStruct<T> ret;
178  ret.values = D;
179  ret.vectors = M;
180  return ret;
181 
182  };
183 
184  //==============================
193  template < typename T >
194  SvdStruct<T> svd(const ublas::matrix<T> & m1)
195  {
196  SvdStruct<T> ret;
197 
198  ublas::matrix<double, ublas::column_major> A(m1.size1(),m1.size2());
199  ublas::matrix<double, ublas::column_major> U(m1.size1(),m1.size2());
200  ublas::matrix<double, ublas::column_major> V(m1.size1(),m1.size2());
201  ublas::vector<double> D(m1.size1());
202 
203  A = m1;
204 
205  lapack::gesdd(A, D, U, V);
206 
207  ret.u = U;
208  ret.v = V;
209  ret.d = D;
210  return ret;
211  };
212 
213  //==============================
220  template < typename T >
221  int rank(const ublas::matrix<T> & m)
222  {
223  SvdStruct<double> svd_value = svd(m);
224  int rank = 0;
225  for(int i=0; i<svd_value.d.size(); i++)
226  {
227  //if singular value is not 'small' (1e-9 is small enough)
228  if(svd_value.d[i] > 1e-9)
229  rank++;
230  }
231  return rank;
232  }
233 
234  //==============================
240  template < typename T >
241  double det(const boost::numeric::ublas::matrix_expression<T>& m) {
242  assert(m().size1() == m().size2());
243 
244  boost::numeric::ublas::permutation_matrix<std::size_t> pivots(m().size1());
245 
246  // create copy
247  boost::numeric::ublas::matrix<typename T::value_type> mLu(m);
248 
249  lu_factorize(mLu, pivots);
250 
251  double detval = 1.0;
252 
253  for (std::size_t i=0; i < pivots.size(); ++i) {
254  if (pivots(i) != i)
255  detval *= -1.0;
256  detval *= mLu(i,i);
257  }
258  return detval;
259  }
260 
261  //==============================
262  //Lag
270  template < typename T >
271  ublas::matrix<T> Lag(const ublas::matrix<T> &data_orig, int n_lags)
272  {
273  ublas::matrix<T> data = data_orig;
274  //check the number of lags
275  if(n_lags < 1)
276  return data;
277 
278  if(n_lags >= ncol(data))
279  {
280  std::cerr << "ERROR: Cannot apply this many lags to data.\n";
281  return data;
282  }
283 
284  //lag the data
285  int rows = nrow(data);
286  int cols = ncol(data);
287 
288  ublas::matrix<T> lagged_data = cppR::createMatrix(0,
289  (n_lags+1) * rows,
290  cols-n_lags);
291  //ublas::matrix<double> d2 = data;
292  for(int i=0; i<= n_lags; i++)
293  {
294  //R Code: laggedData[ ((lag * rows + 1) : ((lag+1) * rows)), ]
295  // <- data[, ((lag+1) : (cols-nLags+lag))];
296  ublas::matrix_slice<ublas::matrix<double> > m1 =
297  submatrix(lagged_data, i * rows,((i+1) * rows) - 1,0,0);
298  ublas::matrix_slice<ublas::matrix<double> > m2 =
299  submatrix(data, 0, 0, i, (cols-n_lags+i) -1);
300  m1 = m2;
301  }
302 
303  return lagged_data;
304  }
305 
306  //==============================
307  //complex matrices
308 
314  template <typename T>
315  ublas::matrix<T>
316  Re(ublas::matrix<std::complex<T> > m)
317  {
318  ublas::matrix<T> ret(m.size1(),m.size2());
319 
320  for(int i=0;i<ret.size1();i++)
321  for(int j=0;j<ret.size2();j++)
322  ret(i,j) = real(m(i,j));
323  return ret;
324  }
325 
331  template <typename T>
332  ublas::matrix<std::complex<T> >
333  Conj(ublas::matrix<std::complex<T> > m)
334  {
335  ublas::matrix<complex<double> > ret(m.size1(),m.size2());
336  for(int i=0;i<ret.size1();i++)
337  for(int j=0;j<ret.size2();j++)
338  ret(i,j) = conj(m(i,j));
339  return ret;
340  }
341 
342 
343  //==============================
344  //FAST FOURIER TRANSFORM
345 
346 #ifdef FFTW3_H
347 
349 
362  ublas::matrix<std::complex<double> >
363  fft(const ublas::matrix<complex<double> > &m, bool inverse=false)
364  {
365  ublas::matrix<complex<double> > ret(m.size1(), m.size2());
366 
367  int rows,cols;
368  rows = m.size1();
369  cols = m.size2();
370 
371  //create in and out matrices
372  fftw_complex *in,*out;
373  in = (fftw_complex*) fftw_malloc(rows * cols * sizeof(fftw_complex));
374  out = (fftw_complex*) fftw_malloc(rows * cols * sizeof(fftw_complex));
375 
376  //create plan
377  fftw_plan p;
378  if(!inverse)
379  {
380  p = fftw_plan_dft_2d(rows,cols,
381  in, out,
382  FFTW_FORWARD, FFTW_ESTIMATE | FFTW_DESTROY_INPUT);
383  }
384  else
385  {
386  p = fftw_plan_dft_2d(rows,cols,
387  in, out,
388  FFTW_BACKWARD, FFTW_ESTIMATE | FFTW_DESTROY_INPUT);
389  }
390  //copy input data into in
391  for(int i=0; i<rows; i++)
392  for(int j=0; j<cols; j++)
393  {
394  in[j+cols*i][0] = real(m(i,j));
395  in[j+cols*i][1] = imag(m(i,j));
396  }
397 
398  //run the fft
399  fftw_execute(p);
400 
401  //copy output data into return matrix
402  for(int i=0; i<rows; i++)
403  for(int j=0; j<cols; j++)
404  {
405  ret(i,j) = complex<double>(out[j+cols*i][0], out[j+cols*i][1]);
406  }
407 
408  //destroy plan
409  fftw_destroy_plan(p);
410 
411  //free arrays
412  fftw_free(in);
413  fftw_free(out);
414 
415  //return matrix
416  return ret;
417  }
418 
419 
429  ublas::matrix<std::complex<double> >
430  fft(const ublas::matrix<double> &m, bool inverse=false)
431  {
432  int rows,cols;
433  rows = m.size1();
434  cols = m.size2();
435 
436  ublas::matrix<complex<double> > ret(rows,cols);
437 
438  for(int i=0; i<rows; i++)
439  for(int j=0; j<cols; j++)
440  {
441  ret(i,j) = complex<double>(m(i,j),0.0);
442  }
443 
444  return fft(ret);
445  }
446 #endif //endif for FFTW
447  //END OF FFT
448 
449  //==============================
450  //OPERATORS
451 
453  template < typename T >
454  ublas::matrix<T>
455  operator+=(const ublas::matrix<T> & m, const T v)
456  {
457  ublas::matrix<T> ret = m;
458  for(int i=0;i<m.size1();i++)
459  for(int j=0;j<m.size2();j++)
460  ret(i,j) += v;
461  return ret;
462  }
463 
465  template < typename T >
466  ublas::matrix<T>
467  operator+(const ublas::matrix<T> & m, const T v)
468  {
469  ublas::matrix<T> ret = m;
470  for(int i=0;i<m.size1();i++)
471  for(int j=0;j<m.size2();j++)
472  ret(i,j) += v;
473  return ret;
474  }
475 
477  template < typename T >
478  ublas::matrix<T>
479  operator-=(const ublas::matrix<T> & m, const T v)
480  {
481  ublas::matrix<T> ret = m;
482  for(int i=0;i<m.size1();i++)
483  for(int j=0;j<m.size2();j++)
484  ret(i,j) -= v;
485  return ret;
486  }
487 
489  template < typename T >
490  ublas::matrix<T>
491  operator-(const ublas::matrix<T> & m, const T v)
492  {
493  ublas::matrix<T> ret = m;
494  for(int i=0;i<m.size1();i++)
495  for(int j=0;j<m.size2();j++)
496  ret(i,j) -= v;
497  return ret;
498  }
499 
500 }//end of namespace
501 
502 
503 //--------------------------------------------------
504 
505 
506 
507 #endif