Nemo  2.4.0b
Simulate forward-in-time genetic evolution in a spatially explicit, individual-based stochastic simulator
tmatrix.h
Go to the documentation of this file.
1 
30 #ifndef _T_MATRIX_H
31 #define _T_MATRIX_H
32 
33 #ifdef HAS_GSL
34 #include <gsl/gsl_vector.h>
35 #include <gsl/gsl_matrix.h>
36 #include <gsl/gsl_permutation.h>
37 #include <gsl/gsl_linalg.h>
38 #endif
39 
40 #include <string.h>
41 #include <sstream>
42 #include <string>
43 #include <assert.h>
44 #include <vector>
45 #include "output.h"
46 
47 using namespace std;
48 
50 class TMatrix {
51 private:
52 
53  unsigned int _rows, _cols, _length;
54 
55  double* _val;
56 
57 public:
58 
59  TMatrix () : _rows(0), _cols(0), _length(0), _val(0) { }
60 
62  TMatrix (const TMatrix& mat) : _rows(0), _cols(0), _length(0), _val(0)
63  {
64  copy(mat);
65  }
66 
68  TMatrix ( unsigned int rows, unsigned int cols ) : _rows(0), _cols(0), _length(0), _val(0)
69  {
70  _length = rows*cols;
71  _val = new double [_length];
72  _rows = rows; _cols = cols;
73  }
74 
75  ~TMatrix () {if(_val != NULL) delete [] _val;}
76 
78  void copy (const TMatrix& mat)
79  {
80  _rows = mat._rows;
81  _cols = mat._cols;
82  _length = _rows*_cols;
83  if(_val) delete [] _val;
84  _val = new double [_length];
85  memcpy(_val,mat._val,_length*sizeof(double));
86  }
87 
90  void copy_recycle (const TMatrix& mat)
91  {
92  if(_rows == mat._rows && _cols == mat._cols)
93  copy(mat);
94  else if(_val == NULL || _length == 0)
95  error("TMatrix::copy_recycle::matrix must be allocated before call\n");
96  else {
97  for(unsigned int i = 0; i < _rows; ++i)
98  for(unsigned int j = 0; j < _cols; ++j)
99  set(i, j, mat._val[ (i%mat._rows)*mat._cols + (j%mat._cols) ]); //mat.get( i%nrow, j%ncol )); _val[i*_cols + j]
100  }
101  }
103  void set (unsigned int i, unsigned int j, double val) {
104  if( i*j < _length)
105  _val[i*_cols + j] = val;
106  else
107  error("TMatrix::set overflow, i*j > length!\n");
108  }
109 
111  void set_row (unsigned int i, double val) {
112  if( i < _rows)
113  for(unsigned int j = 0, stride = i*_cols; j < _cols; ++j)
114  _val[stride + j] = val;
115  else
116  error("TMatrix::set_row overflow, i > num rows!\n");
117  }
118 
120  void set_row (unsigned int i, vector<double> vec) {
121 
122  if(vec.size() < _cols)
123  fatal("TMatrix::set_row copy from vector too small to fit\n");
124 
125  if( i < _rows)
126  for(unsigned int j = 0, stride = i*_cols; j < _cols; ++j)
127  _val[stride + j] = vec[j];
128  else
129  error("TMatrix::set_row overflow, i > num rows!\n");
130  }
131 
133  void set_col (unsigned int i, double val) {
134  if( i < _cols)
135  for(unsigned int j = 0; j < _rows; ++j)
136  _val[j*_cols + i] = val;
137  else
138  error("TMatrix::set_col overflow, i > num cols!\n");
139  }
140 
142  void set_col (unsigned int i, vector<double> vec) {
143 
144  if(vec.size() < _rows)
145  fatal("TMatrix::set_col copy from vector too small to fit\n");
146 
147  if( i < _cols)
148  for(unsigned int j = 0; j < _rows; ++j)
149  _val[j*_cols + i] = vec[j];
150  else
151  error("TMatrix::set_col overflow, i > num cols!\n");
152  }
153 
155  void assign (double val)
156  {
157  for(unsigned int i = 0; i < _length; ++i) _val[i] = val;
158  }
159 
161  void reset (unsigned int rows, unsigned int cols) {
162  _length = rows * cols;
163  if(_length == 0) error("TMatrix::attempt to reset a matrix with size = 0!\n");
164  if(_val != NULL) delete [] _val;
165  _val = new double [_length];
166  memset(_val, 0, _length*sizeof(double));
167  _rows = rows; _cols = cols;
168  }
169 
171  void reset (unsigned int rows, unsigned int cols, double value) {
172  reset(rows, cols);
173  assign(value);
174  }
175 
177  void reset (unsigned int rows, unsigned int cols, const double* array) {
178  reset(rows, cols);
179  memcpy(_val, array, _length * sizeof(double));
180  }
181 
183  void reset ( )
184  {
185  _rows = 0;
186  _cols = 0;
187  _length = 0;
188  if(_val != NULL) delete [] _val;
189  _val = NULL;
190  }
191 
193  double get (unsigned int i, unsigned int j) const {
194  if( !((i+1)*(j+1) > _length) )
195  return _val[i*_cols + j];
196  else
197  fatal("TMatrix::get overflow!\n");
198  return 0;
199  }
201  double* get () const {return _val;}
202  double* getValArray () const {return _val;}
207  unsigned int get_dims (unsigned int* dims) const {
208  if(dims != NULL) { dims[0] = _rows; dims[1] = _cols; }
209  return _length;
210  }
212  unsigned int getNbRows ( ) const {return _rows;}
213  unsigned int nrows () const {return _rows;}
215  unsigned int getNbCols ( ) const {return _cols;}
216  unsigned int ncols () const {return _cols;}
218  unsigned int length ( ) const {return _length;}
224  void getColumnView (unsigned int col, unsigned int n, double* array)
225  {
226  if(col > _cols-1) {
227  error("TMatrix::getColumnView: not that many columns in matrix\n");
228  return;
229  }
230  if(n != _rows) {
231  error("TMatrix::getColumnView: array size not equal to number of rows in matrix\n");
232  return;
233  }
234  for(unsigned int i = 0; i < _rows; ++i)
235  array[i] = _val[i*_cols + col];
236  }
242  void getRowView (unsigned int row, unsigned int n, double* array)
243  {
244  if(row > _rows-1) {
245  error("TMatrix::getRowView: not that many rows in matrix\n");
246  return;
247  }
248  if(n != _cols) {
249  error("TMatrix::getRowView: array size not equal to number of columns in matrix\n");
250  return;
251  }
252  for(unsigned int i = 0, stride = row*_cols; i < _cols; ++i)
253  array[i] = _val[stride + i];
254  }
256  void plus (unsigned int i, unsigned int j, double value)
257  {
258  if( i*j < _length)
259  _val[i*_cols + j] += value;
260  else
261  error("TMatrix::plus overflow!\n");
262  }
264  void matrix_increment (double value)
265  {
266  for(unsigned int i = 0; i < _rows; ++i){
267  for(unsigned int j = 0; j < _cols; ++j){
268  plus(i,j,value);
269  }
270  }
271  }
273  void minus (unsigned int i, unsigned int j, double value)
274  {
275  if( i*j < _length)
276  _val[i*_cols + j] -= value;
277  else
278  error("TMatrix::minus overflow!\n");
279  }
281  void multi (unsigned int i, unsigned int j, double value)
282  {
283  if( i*j < _length)
284  _val[i*_cols + j] *= value;
285  else
286  error("TMatrix::multi overflow!\n");
287  }
289  void divide (unsigned int i, unsigned int j, double value)
290  {
291  if( i*j < _length)
292  _val[i*_cols + j] /= value;
293  else
294  error("TMatrix::divide overflow!\n");
295  }
297  void transpose ()
298  {
299  TMatrix tmp(_cols, _rows);
300 
301  for(unsigned int i = 0; i < _rows; i++)
302  for(unsigned int j = 0; j < _cols; j++)
303  tmp.set(j, i, get(i, j));
304 
305  reset(_cols, _rows, tmp.get());
306  }
308  double colSum (unsigned int col)
309  {
310  assert(col < _cols); //has to be [0, _cols-1]
311  double sum = 0;
312  for (unsigned int i = 0; i < _rows; ++i) {
313  sum += _val[i*_cols + col];
314  }
315  return sum;
316  }
318  double rowSum (unsigned int row)
319  {
320  assert(row < _rows); // has to be [0, _rows-1]
321  double sum = 0;
322  for (unsigned int i = 0; i < _cols; ++i) {
323  sum += _val[row*_cols + i];
324  }
325  return sum;
326  }
327 
328  void show_up()
329  {
330  message("TMatrix dimensions: \nrows = %i, columns = %i, length = %i\n",_rows,_cols, _length);
331  for(unsigned int i = 0; i < _rows; i++) {
332  for(unsigned int j = 0; j < _cols; j++)
333  message("%.3f ",_val[i*_cols + j]);
334  message("\n");
335  }
336  }
338  string to_string ()
339  {
340  ostringstream OUT;
341 
342  OUT << "{";
343  for(unsigned int i = 0; i < _rows; i++) {
344  OUT<<"{";
345  for(unsigned int j = 0; j < _cols-1; j++) {
346  OUT<<_val[i*_cols + j]<<",";
347  }
348  OUT<<_val[i*_cols + _cols-1]<<"}";
349  }
350  OUT << "}";
351 
352  return OUT.str();
353  }
354 
355 #ifdef HAS_GSL
359  void get_gsl_matrix(gsl_matrix* M)
360  {
361  if(M->size1 != _rows)
362  fatal("TMatrix::get_gsl_matrix row size of input matrix doesn't match!\n");
363  if(M->size2 != _cols)
364  fatal("TMatrix::get_gsl_matrix col size of input matrix doesn't match!\n");
365  for(unsigned int i = 0; i < _rows; ++i)
366  for(unsigned int j = 0; j < _cols; ++j)
367  gsl_matrix_set(M,i,j,_val[i*_cols + j]);
368  }
372  void set_from_gsl_matrix (gsl_matrix* M)
373  {
374  if(_val != 0) delete [] _val;
375  _rows = M->size1;
376  _cols = M->size2;
377  _length = _rows*_cols;
378  _val = new double [_length];
379  memcpy(_val, M->data, _length*sizeof(double));
380  }
382  void inverse ()
383  {
384  if(_rows != _cols) {
385  error("TMatrix::inverse matrix must be square to invert!\n");
386  return;
387  }
388 
389  gsl_matrix *S = gsl_matrix_alloc(_rows,_cols);
390  get_gsl_matrix(S);
391  gsl_matrix *LU = gsl_matrix_alloc(_rows,_cols);
392  gsl_permutation *P = gsl_permutation_alloc(_rows);
393  int snum;
394 
395  gsl_matrix_memcpy(LU,S);
396  gsl_linalg_LU_decomp(LU,P,&snum);
397  gsl_linalg_LU_invert(LU,P,S);
398 
399  set_from_gsl_matrix(S);
400 
401  gsl_matrix_free(S);
402  gsl_matrix_free(LU);
403  gsl_permutation_free(P);
404  }
405 #endif
406 };
407 
408 #endif
409 
A class to handle matrix in params, coerces matrix into a vector of same total size.
Definition: tmatrix.h:50
void reset(unsigned int rows, unsigned int cols)
Re-allocate the existing matrix with assigned rows and cols dimensions and all elements to 0.
Definition: tmatrix.h:161
unsigned int getNbRows() const
Gives the number of rows.
Definition: tmatrix.h:212
unsigned int _cols
Definition: tmatrix.h:53
double * get() const
Accessor to the whole array.
Definition: tmatrix.h:201
double * getValArray() const
Definition: tmatrix.h:202
double colSum(unsigned int col)
Sum all elements in a column.
Definition: tmatrix.h:308
void set_col(unsigned int i, vector< double > vec)
Sets element at column i to values stored in vector vec.
Definition: tmatrix.h:142
void set_row(unsigned int i, double val)
Sets all elements at row i to value val.
Definition: tmatrix.h:111
void show_up()
Definition: tmatrix.h:328
unsigned int ncols() const
Definition: tmatrix.h:216
TMatrix()
Definition: tmatrix.h:59
void copy(const TMatrix &mat)
Copy a matrix.
Definition: tmatrix.h:78
unsigned int getNbCols() const
Gives the number of columns.
Definition: tmatrix.h:215
void matrix_increment(double value)
Adds a value to all elements in a matrix.
Definition: tmatrix.h:264
void set(unsigned int i, unsigned int j, double val)
Sets element at row i and column j to value val.
Definition: tmatrix.h:103
double get(unsigned int i, unsigned int j) const
Accessor to element at row i and column j.
Definition: tmatrix.h:193
void getColumnView(unsigned int col, unsigned int n, double *array)
Gives access to a column of the matrix.
Definition: tmatrix.h:224
void set_col(unsigned int i, double val)
Sets element at column i to value val.
Definition: tmatrix.h:133
void assign(double val)
Assigns a value to all element of the matrix.
Definition: tmatrix.h:155
void copy_recycle(const TMatrix &mat)
Copy elements of 'mat', recycling elements of 'mat' if its size is smaller than current matrix.
Definition: tmatrix.h:90
unsigned int _rows
Definition: tmatrix.h:53
void divide(unsigned int i, unsigned int j, double value)
Divide an element of the matrix by a value.
Definition: tmatrix.h:289
void reset()
Reset members to zero state.
Definition: tmatrix.h:183
void transpose()
Transpose the matrix, swaps columns for rows.
Definition: tmatrix.h:297
void reset(unsigned int rows, unsigned int cols, const double *array)
Reset the existing matrix to the new dimensions and copies the array, which has to be of the same tot...
Definition: tmatrix.h:177
double rowSum(unsigned int row)
Sum all elements in a row.
Definition: tmatrix.h:318
void minus(unsigned int i, unsigned int j, double value)
Substracts a value from an element of the matrix.
Definition: tmatrix.h:273
void set_row(unsigned int i, vector< double > vec)
Sets elements at row i to values stored in vector vec.
Definition: tmatrix.h:120
void multi(unsigned int i, unsigned int j, double value)
Multiply an element of the matrix by a value.
Definition: tmatrix.h:281
void plus(unsigned int i, unsigned int j, double value)
Adds a value to an element of the matrix.
Definition: tmatrix.h:256
TMatrix(const TMatrix &mat)
copy constructor.
Definition: tmatrix.h:62
double * _val
Definition: tmatrix.h:55
unsigned int get_dims(unsigned int *dims) const
Accessor to the matrix dimensions.
Definition: tmatrix.h:207
~TMatrix()
Definition: tmatrix.h:75
unsigned int length() const
Returns the number of elements in the matrix.
Definition: tmatrix.h:218
string to_string()
Writes the matrix into a string in Nemo's matrix input format.
Definition: tmatrix.h:338
void getRowView(unsigned int row, unsigned int n, double *array)
Gives access to a row of the matrix.
Definition: tmatrix.h:242
TMatrix(unsigned int rows, unsigned int cols)
Creates an array of doubles of size = rows*cols.
Definition: tmatrix.h:68
void reset(unsigned int rows, unsigned int cols, double value)
Reset the existing matrix to the new dimensions and copies the value to all elements.
Definition: tmatrix.h:171
unsigned int nrows() const
Definition: tmatrix.h:213
void fatal(const char *str,...)
Definition: output.cc:96
int error(const char *str,...)
Definition: output.cc:77
void message(const char *message,...)
Definition: output.cc:40

Generated for Nemo v2.4.0b by  doxygen 1.9.1 -- Nemo is hosted on  Download Nemo

Locations of visitors to this page
Catalogued on GSR