Nemo  2.3.56
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 "output.h"
45
46using namespace std;
47
49class TMatrix {
50private:
51
52 unsigned int _rows, _cols, _length;
53
54 double* _val;
55
56public:
57
58 TMatrix () : _rows(0), _cols(0), _length(0), _val(0) { }
59
61 TMatrix (const TMatrix& mat) : _rows(0), _cols(0), _length(0), _val(0)
62 {
63 copy(mat);
64 }
65
67 TMatrix ( unsigned int rows, unsigned int cols ) : _rows(0), _cols(0), _length(0), _val(0)
68 {
69 _length = rows*cols;
70 _val = new double [_length];
71 _rows = rows; _cols = cols;
72 }
73
74 ~TMatrix () {if(_val != NULL) delete [] _val;}
75
77 void copy (const TMatrix& mat)
78 {
79 _rows = mat._rows;
80 _cols = mat._cols;
82 if(_val) delete [] _val;
83 _val = new double [_length];
84 memcpy(_val,mat._val,_length*sizeof(double));
85 }
86
89 void copy_recycle (const TMatrix& mat)
90 {
91 if(_rows == mat._rows && _cols == mat._cols)
92 copy(mat);
93 else if(_val == NULL || _length == 0)
94 error("TMatrix::copy_recycle::matrix must be allocated before call\n");
95 else {
96 for(unsigned int i = 0; i < _rows; ++i)
97 for(unsigned int j = 0; j < _cols; ++j)
98 set(i, j, mat._val[ (i%mat._rows)*mat._cols + (j%mat._cols) ]); //mat.get( i%nrow, j%ncol )); _val[i*_cols + j]
99 }
100 }
102 void set (unsigned int i, unsigned int j, double val) {
103 if( i*j < _length)
104 _val[i*_cols + j] = val;
105 else
106 error("TMatrix::set overflow!\n");
107 }
108
110 void assign (double val)
111 {
112 for(unsigned int i = 0; i < _length; ++i) _val[i] = val;
113 }
114
116 void reset (unsigned int rows, unsigned int cols) {
117 if(_val != NULL) delete [] _val;
118 _length = rows * cols;
119 _val = new double [_length];
120 memset(_val, 0, _length*sizeof(double));
121 _rows = rows; _cols = cols;
122 }
123
125 void reset (unsigned int rows, unsigned int cols, double value) {
126 reset(rows, cols);
127 assign(value);
128 }
129
131 void reset (unsigned int rows, unsigned int cols, double* array) {
132 reset(rows, cols);
133 memcpy(_val, array, _length * sizeof(double));
134 }
135
137 void reset ( )
138 {
139 _rows = 0;
140 _cols = 0;
141 _length = 0;
142 if(_val != NULL) delete [] _val;
143 _val = NULL;
144 }
145
147 double get (unsigned int i, unsigned int j) {
148 if( !((i+1)*(j+1) > _length) )
149 return _val[i*_cols + j];
150 else
151 fatal("TMatrix::get overflow!\n");
152 return 0;
153 }
155 double* get () const {return _val;}
156 double* getValArray () const {return _val;}
161 unsigned int get_dims (unsigned int* dims) {
162 if(dims != NULL) { dims[0] = _rows; dims[1] = _cols; }
163 return _length;
164 }
166 unsigned int getNbRows ( ) {return _rows;}
167 unsigned int nrows () {return _rows;}
169 unsigned int getNbCols ( ) {return _cols;}
170 unsigned int ncols () {return _cols;}
172 unsigned int length ( ) {return _length;}
178 void getColumnView (unsigned int col, unsigned int n, double* array)
179 {
180 if(col > _cols-1) {
181 error("TMatrix::getColumnView: not that many columns in matrix\n");
182 return;
183 }
184 if(n != _rows) {
185 error("TMatrix::getColumnView: array size not equal to number of rows in matrix\n");
186 return;
187 }
188 for(unsigned int i = 0; i < _rows; ++i)
189 array[i] = _val[i*_cols + col];
190 }
196 void getRowView (unsigned int row, unsigned int n, double* array)
197 {
198 if(row > _rows-1) {
199 error("TMatrix::getRowView: not that many rows in matrix\n");
200 return;
201 }
202 if(n != _cols) {
203 error("TMatrix::getRowView: array size not equal to number of columns in matrix\n");
204 return;
205 }
206 for(unsigned int i = 0, stride = row*_cols; i < _cols; ++i)
207 array[i] = _val[stride + i];
208 }
210 void plus (unsigned int i, unsigned int j, double value)
211 {
212 if( i*j < _length)
213 _val[i*_cols + j] += value;
214 else
215 error("TMatrix::plus overflow!\n");
216 }
218 void matrix_increment (double value)
219 {
220 for(unsigned int i = 0; i < _rows; ++i){
221 for(unsigned int j = 0; j < _cols; ++j){
222 plus(i,j,value);
223 }
224 }
225 }
227 void minus (unsigned int i, unsigned int j, double value)
228 {
229 if( i*j < _length)
230 _val[i*_cols + j] -= value;
231 else
232 error("TMatrix::minus overflow!\n");
233 }
235 void multi (unsigned int i, unsigned int j, double value)
236 {
237 if( i*j < _length)
238 _val[i*_cols + j] *= value;
239 else
240 error("TMatrix::multi overflow!\n");
241 }
243 void divide (unsigned int i, unsigned int j, double value)
244 {
245 if( i*j < _length)
246 _val[i*_cols + j] /= value;
247 else
248 error("TMatrix::divide overflow!\n");
249 }
251 void transpose ()
252 {
253 TMatrix tmp(_cols, _rows);
254
255 for(unsigned int i = 0; i < _rows; i++)
256 for(unsigned int j = 0; j < _cols; j++)
257 tmp.set(j, i, get(i, j));
258
259 reset(_cols, _rows, tmp.get());
260 }
262 double colSum (unsigned int col)
263 {
264 assert(col < _cols); //has to be [0, _cols-1]
265 double sum = 0;
266 for (unsigned int i = 0; i < _rows; ++i) {
267 sum += _val[i*_cols + col];
268 }
269 return sum;
270 }
272 double rowSum (unsigned int row)
273 {
274 assert(row < _rows); // has to be [0, _rows-1]
275 double sum = 0;
276 for (unsigned int i = 0; i < _cols; ++i) {
277 sum += _val[row*_cols + i];
278 }
279 return sum;
280 }
281
282 void show_up()
283 {
284 message("TMatrix dimensions: rows = %i, columns = %i\n",_rows,_cols);
285 for(unsigned int i = 0; i < _rows; i++) {
286 for(unsigned int j = 0; j < _cols; j++)
287 message("%.3f ",_val[i*_cols + j]);
288 message("\n");
289 }
290 }
292 string to_string ()
293 {
294 ostringstream OUT;
295
296 OUT << "{";
297 for(unsigned int i = 0; i < _rows; i++) {
298 OUT<<"{";
299 for(unsigned int j = 0; j < _cols-1; j++) {
300 OUT<<_val[i*_cols + j]<<",";
301 }
302 OUT<<_val[i*_cols + _cols-1]<<"}";
303 }
304 OUT << "}";
305
306 return OUT.str();
307 }
308
309#ifdef HAS_GSL
313 void get_gsl_matrix(gsl_matrix* M)
314 {
315 if(M->size1 != _rows)
316 fatal("TMatrix::get_gsl_matrix row size of input matrix doesn't match!\n");
317 if(M->size2 != _cols)
318 fatal("TMatrix::get_gsl_matrix col size of input matrix doesn't match!\n");
319 for(unsigned int i = 0; i < _rows; ++i)
320 for(unsigned int j = 0; j < _cols; ++j)
321 gsl_matrix_set(M,i,j,_val[i*_cols + j]);
322 }
326 void set_from_gsl_matrix (gsl_matrix* M)
327 {
328 if(_val != 0) delete [] _val;
329 _rows = M->size1;
330 _cols = M->size2;
332 _val = new double [_length];
333 memcpy(_val, M->data, _length*sizeof(double));
334 }
336 void inverse ()
337 {
338 if(_rows != _cols) {
339 error("TMatrix::inverse matrix must be square to invert!\n");
340 return;
341 }
342
343 gsl_matrix *S = gsl_matrix_alloc(_rows,_cols);
344 get_gsl_matrix(S);
345 gsl_matrix *LU = gsl_matrix_alloc(_rows,_cols);
346 gsl_permutation *P = gsl_permutation_alloc(_rows);
347 int snum;
348
349 gsl_matrix_memcpy(LU,S);
350 gsl_linalg_LU_decomp(LU,P,&snum);
351 gsl_linalg_LU_invert(LU,P,S);
352
353 set_from_gsl_matrix(S);
354
355 gsl_matrix_free(S);
356 gsl_matrix_free(LU);
357 gsl_permutation_free(P);
358 }
359#endif
360};
361
362#endif
363
A class to handle matrix in params, coerces matrix into a vector of same total size.
Definition: tmatrix.h:49
void reset(unsigned int rows, unsigned int cols)
Re-allocate the existing matrix with assigned rows and cols dimensions.
Definition: tmatrix.h:116
unsigned int _cols
Definition: tmatrix.h:52
unsigned int length()
Returns the number of elements in the matrix.
Definition: tmatrix.h:172
double * getValArray() const
Definition: tmatrix.h:156
double colSum(unsigned int col)
Sum all elements in a column.
Definition: tmatrix.h:262
void show_up()
Definition: tmatrix.h:282
TMatrix()
Definition: tmatrix.h:58
unsigned int ncols()
Definition: tmatrix.h:170
void copy(const TMatrix &mat)
Copy a matrix.
Definition: tmatrix.h:77
void matrix_increment(double value)
Adds a value to all elements in a matrix.
Definition: tmatrix.h:218
void set(unsigned int i, unsigned int j, double val)
Sets element at row i and column j to value val.
Definition: tmatrix.h:102
unsigned int get_dims(unsigned int *dims)
Accessor to the matrix dimensions.
Definition: tmatrix.h:161
void getColumnView(unsigned int col, unsigned int n, double *array)
Gives access to a column of the matrix.
Definition: tmatrix.h:178
void reset(unsigned int rows, unsigned int cols, double *array)
Reset the existing matrix to the new dimensions and copies the array.
Definition: tmatrix.h:131
void assign(double val)
Assigns a value to all element of the matrix.
Definition: tmatrix.h:110
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:89
unsigned int _rows
Definition: tmatrix.h:52
void divide(unsigned int i, unsigned int j, double value)
Divide an element of the matrix by a value.
Definition: tmatrix.h:243
void reset()
Reset members to zero state.
Definition: tmatrix.h:137
void transpose()
Transpose the matrix, swaps columns for rows.
Definition: tmatrix.h:251
double rowSum(unsigned int row)
Sum all elements in a row.
Definition: tmatrix.h:272
void minus(unsigned int i, unsigned int j, double value)
Substracts a value from an element of the matrix.
Definition: tmatrix.h:227
double * get() const
Accessor to the whole array.
Definition: tmatrix.h:155
double get(unsigned int i, unsigned int j)
Accessor to element at row i and column j.
Definition: tmatrix.h:147
unsigned int getNbRows()
Gives the number of rows.
Definition: tmatrix.h:166
unsigned int nrows()
Definition: tmatrix.h:167
void multi(unsigned int i, unsigned int j, double value)
Multiply an element of the matrix by a value.
Definition: tmatrix.h:235
void plus(unsigned int i, unsigned int j, double value)
Adds a value to an element of the matrix.
Definition: tmatrix.h:210
TMatrix(const TMatrix &mat)
copy constructor.
Definition: tmatrix.h:61
double * _val
Definition: tmatrix.h:54
unsigned int _length
Definition: tmatrix.h:52
~TMatrix()
Definition: tmatrix.h:74
string to_string()
Writes the matrix into a string in Nemo's matrix input format.
Definition: tmatrix.h:292
void getRowView(unsigned int row, unsigned int n, double *array)
Gives access to a row of the matrix.
Definition: tmatrix.h:196
TMatrix(unsigned int rows, unsigned int cols)
Creates an array of doubles of size = rows*cols.
Definition: tmatrix.h:67
unsigned int getNbCols()
Gives the number of columns.
Definition: tmatrix.h:169
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:125
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.3.56 by  doxygen 1.9.0 -- Nemo is hosted on  Download Nemo

Locations of visitors to this page
Catalogued on GSR