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