1 | // $Id: matrix.cc 916 2007-09-30 00:50:10Z peter $ |
---|
2 | |
---|
3 | /* |
---|
4 | Copyright (C) 2003 Daniel Dalevi, Peter Johansson |
---|
5 | Copyright (C) 2004 Jari Häkkinen, Peter Johansson |
---|
6 | Copyright (C) 2005, 2006 Jari Häkkinen, Markus Ringnér, Peter Johansson |
---|
7 | Copyright (C) 2007 Jari Häkkinen, Peter Johansson |
---|
8 | |
---|
9 | This file is part of the yat library, http://trac.thep.lu.se/trac/yat |
---|
10 | |
---|
11 | The yat library is free software; you can redistribute it and/or |
---|
12 | modify it under the terms of the GNU General Public License as |
---|
13 | published by the Free Software Foundation; either version 2 of the |
---|
14 | License, or (at your option) any later version. |
---|
15 | |
---|
16 | The yat library is distributed in the hope that it will be useful, |
---|
17 | but WITHOUT ANY WARRANTY; without even the implied warranty of |
---|
18 | MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU |
---|
19 | General Public License for more details. |
---|
20 | |
---|
21 | You should have received a copy of the GNU General Public License |
---|
22 | along with this program; if not, write to the Free Software |
---|
23 | Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA |
---|
24 | 02111-1307, USA. |
---|
25 | */ |
---|
26 | |
---|
27 | #include "matrix.h" |
---|
28 | #include "vector.h" |
---|
29 | #include "utility.h" |
---|
30 | |
---|
31 | #include <cassert> |
---|
32 | #include <cmath> |
---|
33 | #include <sstream> |
---|
34 | #include <vector> |
---|
35 | |
---|
36 | #include <gsl/gsl_blas.h> |
---|
37 | |
---|
38 | namespace theplu { |
---|
39 | namespace yat { |
---|
40 | namespace utility { |
---|
41 | |
---|
42 | |
---|
43 | matrix::matrix(void) |
---|
44 | : blas_result_(NULL), m_(NULL), view_(NULL), view_const_(NULL), |
---|
45 | proxy_m_(NULL) |
---|
46 | { |
---|
47 | } |
---|
48 | |
---|
49 | |
---|
50 | matrix::matrix(const size_t& r, const size_t& c, double init_value) |
---|
51 | : blas_result_(NULL), m_(gsl_matrix_alloc(r,c)), view_(NULL), |
---|
52 | view_const_(NULL), proxy_m_(m_) |
---|
53 | { |
---|
54 | if (!m_) |
---|
55 | throw utility::GSL_error("matrix::matrix failed to allocate memory"); |
---|
56 | all(init_value); |
---|
57 | } |
---|
58 | |
---|
59 | |
---|
60 | matrix::matrix(const matrix& o) |
---|
61 | : blas_result_(NULL), m_(o.create_gsl_matrix_copy()), view_(NULL), |
---|
62 | view_const_(NULL), proxy_m_(m_) |
---|
63 | { |
---|
64 | } |
---|
65 | |
---|
66 | |
---|
67 | matrix::matrix(matrix& m, size_t offset_row, size_t offset_column, |
---|
68 | size_t n_row, size_t n_column) |
---|
69 | : blas_result_(NULL), view_const_(NULL) |
---|
70 | { |
---|
71 | view_ = new gsl_matrix_view(gsl_matrix_submatrix(m.m_, |
---|
72 | offset_row,offset_column, |
---|
73 | n_row,n_column)); |
---|
74 | if (!view_) |
---|
75 | throw utility::GSL_error("matrix::matrix failed to setup view"); |
---|
76 | proxy_m_ = m_ = &(view_->matrix); |
---|
77 | } |
---|
78 | |
---|
79 | |
---|
80 | // Constructor that gets data from istream |
---|
81 | matrix::matrix(std::istream& is, char sep) |
---|
82 | throw (utility::IO_error,std::exception) |
---|
83 | : blas_result_(NULL), view_(NULL), view_const_(NULL) |
---|
84 | { |
---|
85 | // read the data file and store in stl vectors (dynamically |
---|
86 | // expandable) |
---|
87 | std::vector<std::vector<double> > data_matrix; |
---|
88 | u_int nof_columns=0; |
---|
89 | u_int nof_rows = 0; |
---|
90 | std::string line; |
---|
91 | while(getline(is, line, '\n')){ |
---|
92 | // Ignoring empty lines |
---|
93 | if (!line.size()) { |
---|
94 | continue; |
---|
95 | } |
---|
96 | nof_rows++; |
---|
97 | std::vector<double> v; |
---|
98 | std::string element; |
---|
99 | std::stringstream ss(line); |
---|
100 | |
---|
101 | bool ok=true; |
---|
102 | while(ok) { |
---|
103 | if(sep=='\0') |
---|
104 | ok=(ss>>element); |
---|
105 | else |
---|
106 | ok=getline(ss, element, sep); |
---|
107 | if(!ok) |
---|
108 | break; |
---|
109 | |
---|
110 | if(utility::is_double(element)) { |
---|
111 | v.push_back(atof(element.c_str())); |
---|
112 | } |
---|
113 | else if (!element.size() || utility::is_nan(element)) { |
---|
114 | v.push_back(std::numeric_limits<double>::quiet_NaN()); |
---|
115 | } |
---|
116 | else { |
---|
117 | std::stringstream ss("Warning: '"); |
---|
118 | ss << element << "' is not accepted as a matrix element."; |
---|
119 | throw IO_error(ss.str()); |
---|
120 | } |
---|
121 | } |
---|
122 | if(sep!='\0' && line[line.size()-1]==sep) // add NaN for final separator |
---|
123 | v.push_back(std::numeric_limits<double>::quiet_NaN()); |
---|
124 | if (!nof_columns) |
---|
125 | nof_columns=v.size(); |
---|
126 | else if (v.size()!=nof_columns) { |
---|
127 | std::ostringstream s; |
---|
128 | s << "matrix::matrix(std::istream&, char) data file error: " |
---|
129 | << "line " << nof_rows << " has " << v.size() |
---|
130 | << " columns; expected " << nof_columns << " columns."; |
---|
131 | throw utility::IO_error(s.str()); |
---|
132 | } |
---|
133 | data_matrix.push_back(v); |
---|
134 | } |
---|
135 | |
---|
136 | // manipulate the state of the stream to be good |
---|
137 | is.clear(std::ios::goodbit); |
---|
138 | // convert the data to a gsl matrix |
---|
139 | proxy_m_ = m_ = gsl_matrix_alloc ( nof_rows, nof_columns ); |
---|
140 | if (!m_) |
---|
141 | throw utility::GSL_error("matrix::matrix failed to allocate memory"); |
---|
142 | |
---|
143 | // if gsl error handler disabled, out of bounds index will not |
---|
144 | // abort the program. |
---|
145 | for(u_int i=0;i<nof_rows;i++) |
---|
146 | for(u_int j=0;j<nof_columns;j++) |
---|
147 | gsl_matrix_set( m_, i, j, data_matrix[i][j] ); |
---|
148 | } |
---|
149 | |
---|
150 | |
---|
151 | matrix::~matrix(void) |
---|
152 | { |
---|
153 | delete_allocated_memory(); |
---|
154 | if (blas_result_) |
---|
155 | gsl_matrix_free(blas_result_); |
---|
156 | } |
---|
157 | |
---|
158 | |
---|
159 | const matrix& matrix::clone(const matrix& other) |
---|
160 | { |
---|
161 | if (this!=&other) { |
---|
162 | |
---|
163 | delete_allocated_memory(); |
---|
164 | |
---|
165 | if (other.view_) { |
---|
166 | view_ = new gsl_matrix_view(*other.view_); |
---|
167 | proxy_m_ = m_ = &(view_->matrix); |
---|
168 | } |
---|
169 | else if (other.view_const_) { |
---|
170 | view_const_ = new gsl_matrix_const_view(*other.view_const_); |
---|
171 | proxy_m_ = &(view_const_->matrix); |
---|
172 | } |
---|
173 | else if (other.m_) |
---|
174 | proxy_m_ = m_ = other.create_gsl_matrix_copy(); |
---|
175 | |
---|
176 | // no need to delete blas_result_ if the number of rows fit, it |
---|
177 | // may be useful later. |
---|
178 | if (blas_result_ && (blas_result_->size1!=rows())) { |
---|
179 | gsl_matrix_free(blas_result_); |
---|
180 | blas_result_=NULL; |
---|
181 | } |
---|
182 | } |
---|
183 | return *this; |
---|
184 | } |
---|
185 | |
---|
186 | |
---|
187 | size_t matrix::columns(void) const |
---|
188 | { |
---|
189 | if (!proxy_m_) |
---|
190 | return 0; |
---|
191 | return proxy_m_->size2; |
---|
192 | } |
---|
193 | |
---|
194 | |
---|
195 | gsl_matrix* matrix::create_gsl_matrix_copy(void) const |
---|
196 | { |
---|
197 | gsl_matrix* m = gsl_matrix_alloc(rows(),columns()); |
---|
198 | if (!m) |
---|
199 | throw utility::GSL_error("matrix::create_gsl_matrix_copy failed to allocate memory"); |
---|
200 | if (gsl_matrix_memcpy(m,proxy_m_)) |
---|
201 | throw utility::GSL_error("matrix::create_gsl_matrix_copy dimension mis-match"); |
---|
202 | return m; |
---|
203 | } |
---|
204 | |
---|
205 | |
---|
206 | void matrix::delete_allocated_memory(void) |
---|
207 | { |
---|
208 | if (view_) |
---|
209 | delete view_; |
---|
210 | else if (view_const_) |
---|
211 | delete view_const_; |
---|
212 | else if (m_) |
---|
213 | gsl_matrix_free(m_); |
---|
214 | blas_result_=NULL; |
---|
215 | proxy_m_=m_=NULL; |
---|
216 | } |
---|
217 | |
---|
218 | |
---|
219 | void matrix::div(const matrix& other) |
---|
220 | { |
---|
221 | assert(m_); |
---|
222 | int status=gsl_matrix_div_elements(m_, other.gsl_matrix_p()); |
---|
223 | if (status) |
---|
224 | throw utility::GSL_error(std::string("matrix::div_elements",status)); |
---|
225 | } |
---|
226 | |
---|
227 | |
---|
228 | bool matrix::equal(const matrix& other, const double d) const |
---|
229 | { |
---|
230 | if (this==&other) |
---|
231 | return true; |
---|
232 | if (columns()!=other.columns() || rows()!=other.rows()) |
---|
233 | return false; |
---|
234 | for (size_t i=0; i<rows(); i++) |
---|
235 | for (size_t j=0; j<columns(); j++) |
---|
236 | // The two last condition checks are needed for NaN detection |
---|
237 | if (fabs( (*this)(i,j)-other(i,j) ) > d || |
---|
238 | (*this)(i,j)!=(*this)(i,j) || other(i,j)!=other(i,j)) |
---|
239 | return false; |
---|
240 | return true; |
---|
241 | } |
---|
242 | |
---|
243 | |
---|
244 | const gsl_matrix* matrix::gsl_matrix_p(void) const |
---|
245 | { |
---|
246 | return proxy_m_; |
---|
247 | } |
---|
248 | |
---|
249 | |
---|
250 | gsl_matrix* matrix::gsl_matrix_p(void) |
---|
251 | { |
---|
252 | return m_; |
---|
253 | } |
---|
254 | |
---|
255 | |
---|
256 | bool matrix::isview(void) const |
---|
257 | { |
---|
258 | return view_ || view_const_; |
---|
259 | } |
---|
260 | |
---|
261 | |
---|
262 | void matrix::mul(const matrix& other) |
---|
263 | { |
---|
264 | assert(m_); |
---|
265 | int status=gsl_matrix_mul_elements(m_, other.gsl_matrix_p()); |
---|
266 | if (status) |
---|
267 | throw utility::GSL_error(std::string("matrix::mul_elements",status)); |
---|
268 | } |
---|
269 | |
---|
270 | |
---|
271 | void matrix::resize(size_t r, size_t c, double init_value) |
---|
272 | { |
---|
273 | delete_allocated_memory(); |
---|
274 | |
---|
275 | proxy_m_ = m_ = gsl_matrix_alloc(r,c); |
---|
276 | if (!m_) |
---|
277 | throw utility::GSL_error("matrix::matrix failed to allocate memory"); |
---|
278 | all(init_value); |
---|
279 | |
---|
280 | // no need to delete blas_result_ if the number of rows fit, it |
---|
281 | // may be useful later. |
---|
282 | if (blas_result_ && (blas_result_->size1!=rows())) { |
---|
283 | gsl_matrix_free(blas_result_); |
---|
284 | blas_result_=NULL; |
---|
285 | } |
---|
286 | } |
---|
287 | |
---|
288 | |
---|
289 | size_t matrix::rows(void) const |
---|
290 | { |
---|
291 | if (!proxy_m_) |
---|
292 | return 0; |
---|
293 | return proxy_m_->size1; |
---|
294 | } |
---|
295 | |
---|
296 | |
---|
297 | void matrix::all(const double value) |
---|
298 | { |
---|
299 | assert(m_); |
---|
300 | gsl_matrix_set_all(m_, value); |
---|
301 | } |
---|
302 | |
---|
303 | |
---|
304 | void matrix::column(const size_t column, const vector& vec) |
---|
305 | { |
---|
306 | assert(m_); |
---|
307 | int status=gsl_matrix_set_col(m_, column, vec.gsl_vector_p()); |
---|
308 | if (status) |
---|
309 | throw utility::GSL_error(std::string("matrix::set_column",status)); |
---|
310 | } |
---|
311 | |
---|
312 | |
---|
313 | void matrix::row(const size_t row, const vector& vec) |
---|
314 | { |
---|
315 | assert(m_); |
---|
316 | int status=gsl_matrix_set_row(m_, row, vec.gsl_vector_p()); |
---|
317 | if (status) |
---|
318 | throw utility::GSL_error(std::string("matrix::set_row",status)); |
---|
319 | } |
---|
320 | |
---|
321 | |
---|
322 | void matrix::swap_columns(const size_t i, const size_t j) |
---|
323 | { |
---|
324 | assert(m_); |
---|
325 | int status=gsl_matrix_swap_columns(m_, i, j); |
---|
326 | if (status) |
---|
327 | throw utility::GSL_error(std::string("matrix::swap_columns",status)); |
---|
328 | } |
---|
329 | |
---|
330 | |
---|
331 | void matrix::swap_rowcol(const size_t i, const size_t j) |
---|
332 | { |
---|
333 | assert(m_); |
---|
334 | int status=gsl_matrix_swap_rowcol(m_, i, j); |
---|
335 | if (status) |
---|
336 | throw utility::GSL_error(std::string("matrix::swap_rowcol",status)); |
---|
337 | } |
---|
338 | |
---|
339 | |
---|
340 | void matrix::swap_rows(const size_t i, const size_t j) |
---|
341 | { |
---|
342 | assert(m_); |
---|
343 | int status=gsl_matrix_swap_rows(m_, i, j); |
---|
344 | if (status) |
---|
345 | throw utility::GSL_error(std::string("matrix::swap_rows",status)); |
---|
346 | } |
---|
347 | |
---|
348 | |
---|
349 | void matrix::transpose(void) |
---|
350 | { |
---|
351 | assert(m_); |
---|
352 | if (columns()==rows()) |
---|
353 | gsl_matrix_transpose(m_); // this never fails |
---|
354 | else { |
---|
355 | gsl_matrix* transposed = gsl_matrix_alloc(columns(),rows()); |
---|
356 | if (!transposed) |
---|
357 | throw utility::GSL_error("matrix::transpose failed to allocate memory"); |
---|
358 | // next line never fails if allocation above succeeded. |
---|
359 | gsl_matrix_transpose_memcpy(transposed,m_); |
---|
360 | gsl_matrix_free(m_); |
---|
361 | proxy_m_ = m_ = transposed; |
---|
362 | if (blas_result_) { |
---|
363 | gsl_matrix_free(blas_result_); |
---|
364 | blas_result_=NULL; |
---|
365 | } |
---|
366 | } |
---|
367 | } |
---|
368 | |
---|
369 | |
---|
370 | double& matrix::operator()(size_t row, size_t column) |
---|
371 | { |
---|
372 | assert(m_); |
---|
373 | assert(row<rows()); |
---|
374 | assert(column<columns()); |
---|
375 | double* d=gsl_matrix_ptr(m_, row, column); |
---|
376 | if (!d) |
---|
377 | throw utility::GSL_error("matrix::operator()",GSL_EINVAL); |
---|
378 | return *d; |
---|
379 | } |
---|
380 | |
---|
381 | |
---|
382 | const double& matrix::operator()(size_t row, size_t column) const |
---|
383 | { |
---|
384 | assert(row<rows()); |
---|
385 | assert(column<columns()); |
---|
386 | const double* d=gsl_matrix_const_ptr(proxy_m_, row, column); |
---|
387 | if (!d) |
---|
388 | throw utility::GSL_error("matrix::operator()",GSL_EINVAL); |
---|
389 | return *d; |
---|
390 | } |
---|
391 | |
---|
392 | |
---|
393 | bool matrix::operator==(const matrix& other) const |
---|
394 | { |
---|
395 | return equal(other); |
---|
396 | } |
---|
397 | |
---|
398 | |
---|
399 | bool matrix::operator!=(const matrix& other) const |
---|
400 | { |
---|
401 | return !equal(other); |
---|
402 | } |
---|
403 | |
---|
404 | |
---|
405 | const matrix& matrix::operator=( const matrix& other ) |
---|
406 | { |
---|
407 | assert(m_); |
---|
408 | if (this!=&other) |
---|
409 | if (gsl_matrix_memcpy(m_, other.gsl_matrix_p())) |
---|
410 | throw utility::GSL_error("matrix::create_gsl_matrix_copy dimension mis-match"); |
---|
411 | return *this; |
---|
412 | } |
---|
413 | |
---|
414 | |
---|
415 | const matrix& matrix::operator+=(const matrix& other) |
---|
416 | { |
---|
417 | assert(m_); |
---|
418 | int status=gsl_matrix_add(m_, other.proxy_m_); |
---|
419 | if (status) |
---|
420 | throw utility::GSL_error(std::string("matrix::operator+=", status)); |
---|
421 | return *this; |
---|
422 | } |
---|
423 | |
---|
424 | |
---|
425 | const matrix& matrix::operator+=(const double d) |
---|
426 | { |
---|
427 | assert(m_); |
---|
428 | gsl_matrix_add_constant(m_, d); |
---|
429 | return *this; |
---|
430 | } |
---|
431 | |
---|
432 | |
---|
433 | const matrix& matrix::operator-=(const matrix& other) |
---|
434 | { |
---|
435 | assert(m_); |
---|
436 | int status=gsl_matrix_sub(m_, other.proxy_m_); |
---|
437 | if (status) |
---|
438 | throw utility::GSL_error(std::string("matrix::operator-=", status)); |
---|
439 | return *this; |
---|
440 | } |
---|
441 | |
---|
442 | |
---|
443 | const matrix& matrix::operator-=(const double d) |
---|
444 | { |
---|
445 | assert(m_); |
---|
446 | gsl_matrix_add_constant(m_, -d); |
---|
447 | return *this; |
---|
448 | } |
---|
449 | |
---|
450 | |
---|
451 | const matrix& matrix::operator*=(const matrix& other) |
---|
452 | { |
---|
453 | assert(m_); |
---|
454 | if ( blas_result_ && ((blas_result_->size1!=rows()) || |
---|
455 | (blas_result_->size2!=other.columns())) ) { |
---|
456 | gsl_matrix_free(blas_result_); |
---|
457 | blas_result_=NULL; |
---|
458 | } |
---|
459 | if (!blas_result_) { |
---|
460 | blas_result_ = gsl_matrix_alloc(rows(),other.columns()); |
---|
461 | if (!blas_result_) |
---|
462 | throw utility::GSL_error("matrix::operator*= failed to allocate memory"); |
---|
463 | } |
---|
464 | gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, m_, other.proxy_m_, 0.0, |
---|
465 | blas_result_); |
---|
466 | gsl_matrix* tmp=m_; |
---|
467 | proxy_m_ = m_ = blas_result_; |
---|
468 | blas_result_=tmp; |
---|
469 | return *this; |
---|
470 | } |
---|
471 | |
---|
472 | |
---|
473 | const matrix& matrix::operator*=(const double d) |
---|
474 | { |
---|
475 | assert(m_); |
---|
476 | gsl_matrix_scale(m_, d); |
---|
477 | return *this; |
---|
478 | } |
---|
479 | |
---|
480 | |
---|
481 | bool isnull(const matrix& other) |
---|
482 | { |
---|
483 | return gsl_matrix_isnull(other.gsl_matrix_p()); |
---|
484 | } |
---|
485 | |
---|
486 | |
---|
487 | double max(const matrix& other) |
---|
488 | { |
---|
489 | return gsl_matrix_max(other.gsl_matrix_p()); |
---|
490 | } |
---|
491 | |
---|
492 | |
---|
493 | double min(const matrix& other) |
---|
494 | { |
---|
495 | return gsl_matrix_min(other.gsl_matrix_p()); |
---|
496 | } |
---|
497 | |
---|
498 | |
---|
499 | void minmax_index(const matrix& other, |
---|
500 | std::pair<size_t,size_t>& min, std::pair<size_t,size_t>& max) |
---|
501 | { |
---|
502 | gsl_matrix_minmax_index(other.gsl_matrix_p(), &min.first, &min.second, |
---|
503 | &max.first, &max.second); |
---|
504 | } |
---|
505 | |
---|
506 | |
---|
507 | bool nan(const matrix& templat, matrix& flag) |
---|
508 | { |
---|
509 | size_t rows=templat.rows(); |
---|
510 | size_t columns=templat.columns(); |
---|
511 | if (rows!=flag.rows() && columns!=flag.columns()) |
---|
512 | flag.clone(matrix(rows,columns,1.0)); |
---|
513 | else |
---|
514 | flag.all(1.0); |
---|
515 | bool nan=false; |
---|
516 | for (size_t i=0; i<rows; i++) |
---|
517 | for (size_t j=0; j<columns; j++) |
---|
518 | if (std::isnan(templat(i,j))) { |
---|
519 | flag(i,j)=0; |
---|
520 | nan=true; |
---|
521 | } |
---|
522 | return nan; |
---|
523 | } |
---|
524 | |
---|
525 | |
---|
526 | void swap(matrix& a, matrix& b) |
---|
527 | { |
---|
528 | assert(a.gsl_matrix_p()); assert(b.gsl_matrix_p()); |
---|
529 | int status=gsl_matrix_swap(a.gsl_matrix_p(), b.gsl_matrix_p()); |
---|
530 | if (status) |
---|
531 | throw utility::GSL_error(std::string("swap(matrix&,matrix&)",status)); |
---|
532 | } |
---|
533 | |
---|
534 | |
---|
535 | std::ostream& operator<<(std::ostream& s, const matrix& m) |
---|
536 | { |
---|
537 | s.setf(std::ios::dec); |
---|
538 | s.precision(12); |
---|
539 | for(size_t i=0, j=0; i<m.rows(); i++) |
---|
540 | for (j=0; j<m.columns(); j++) { |
---|
541 | s << m(i,j); |
---|
542 | if (j<m.columns()-1) |
---|
543 | s << s.fill(); |
---|
544 | else if (i<m.rows()-1) |
---|
545 | s << "\n"; |
---|
546 | } |
---|
547 | return s; |
---|
548 | } |
---|
549 | |
---|
550 | |
---|
551 | vector operator*(const matrix& m, const vector& v) |
---|
552 | { |
---|
553 | utility::vector res(m.rows()); |
---|
554 | for (size_t i=0; i<res.size(); ++i) |
---|
555 | res(i) = vector(m,i) * v; |
---|
556 | return res; |
---|
557 | } |
---|
558 | |
---|
559 | |
---|
560 | vector operator*(const vector& v, const matrix& m) |
---|
561 | { |
---|
562 | utility::vector res(m.columns()); |
---|
563 | for (size_t i=0; i<res.size(); ++i) |
---|
564 | res(i) = v * vector(m,i,false); |
---|
565 | return res; |
---|
566 | } |
---|
567 | |
---|
568 | }}} // of namespace utility, yat and thep |
---|