source: trunk/yat/utility/matrix.cc @ 754

Last change on this file since 754 was 754, checked in by Jari Häkkinen, 16 years ago

Addresses #2 and #65. Continued adding GSL_error exceptions, and added doxygen briefs.

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 10.1 KB
Line 
1// $Id: matrix.cc 754 2007-02-17 22:33:44Z jari $
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 Jari Häkkinen
8
9  This file is part of the yat library, http://lev.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 <cmath>
32#include <sstream>
33#include <vector>
34
35#include <gsl/gsl_blas.h>
36
37namespace theplu {
38namespace yat {
39namespace utility {
40
41
42  matrix::matrix(void)
43    : m_(NULL), view_(NULL)
44  {
45  }
46
47
48  matrix::matrix(const size_t& r, const size_t& c, double init_value)
49    : m_(gsl_matrix_alloc(r,c)), view_(NULL)
50  {
51    if (!m_)
52      throw utility::GSL_error("matrix::matrix failed to allocate memory");
53    set_all(init_value);
54  }
55
56
57  matrix::matrix(const matrix& o)
58    : m_(o.create_gsl_matrix_copy()), view_(NULL)
59  {
60  }
61
62
63  matrix::matrix(matrix& m, size_t offset_row, size_t offset_column,
64                 size_t n_row, size_t n_column)
65  {
66    view_ = new gsl_matrix_view(gsl_matrix_submatrix(m.m_,
67                                                     offset_row,offset_column,
68                                                     n_row,n_column));
69    if (!view_)
70      throw utility::GSL_error("matrix::matrix failed to setup view");
71    m_ = &(view_->matrix);
72  }
73
74
75  // Constructor that gets data from istream
76  matrix::matrix(std::istream& is, char sep) 
77    throw (utility::IO_error,std::exception)
78    : view_(NULL)
79  {
80    // Markus to Jari, somewhere we should check that quiet_NaNs are supported
81    // std::numeric_limits<double>::has_quiet_NaN has to be true.
82    // Also in vector
83
84    // read the data file and store in stl vectors (dynamically
85    // expandable)
86    std::vector<std::vector<double> > data_matrix;
87    u_int nof_columns=0;
88    u_int nof_rows = 0;
89    std::string line;
90    while(getline(is, line, '\n')){
91      // Ignoring empty lines
92      if (!line.size()) {
93        continue;
94      }
95      nof_rows++;
96      std::vector<double> v;
97      std::string element;
98      std::stringstream ss(line);
99     
100      bool ok=true;
101      while(ok) {
102        if(sep=='\0')
103          ok=(ss>>element);
104        else
105          ok=getline(ss, element, sep);
106        if(!ok)
107          break;
108       
109        if(utility::is_double(element)) {
110          v.push_back(atof(element.c_str()));
111        }
112        else if (!element.size() || utility::is_nan(element)) {
113          v.push_back(std::numeric_limits<double>::quiet_NaN());
114        }
115        else {
116          // Jari, this should be communicated with as an exception.
117          std::cerr << "Warning: '" << element
118                    << "' is not accepted as a matrix element." << std::endl;
119        }
120      }           
121      if(sep!='\0' && line[line.size()-1]==sep) // add NaN for final separator
122          v.push_back(std::numeric_limits<double>::quiet_NaN());
123      if (!nof_columns)
124        nof_columns=v.size();
125      else if (v.size()!=nof_columns) {
126        std::ostringstream s;
127        s << "matrix::matrix(std::istream&, char) data file error: "
128          << "line " << nof_rows << " has " << v.size()
129          << " columns; expected " << nof_columns << " columns.";
130        throw utility::IO_error(s.str());
131      }
132      data_matrix.push_back(v);
133    }
134
135    // manipulate the state of the stream to be good
136    is.clear(std::ios::goodbit);
137    // convert the data to a gsl matrix
138    m_ = gsl_matrix_alloc ( nof_rows, nof_columns );
139    if (!m_)
140      throw utility::GSL_error("matrix::matrix failed to allocate memory");
141
142    // if gsl error handler disabled, out of bounds index will not
143    // abort the program.
144    for(u_int i=0;i<nof_rows;i++)
145      for(u_int j=0;j<nof_columns;j++)
146        gsl_matrix_set( m_, i, j, data_matrix[i][j] );
147  }
148
149
150  matrix::~matrix(void)
151  {
152    if (view_)
153      delete view_;
154    else if (m_)
155      gsl_matrix_free(m_);
156    m_=NULL;
157  }
158
159
160  void matrix::add(const matrix& b)
161  {
162    int status=gsl_matrix_add(m_, b.m_);
163    if (status)
164      throw utility::GSL_error(std::string("matrix::add",status));
165  }
166
167
168  void matrix::add(double d)
169  {
170    gsl_matrix_add_constant(m_, d);
171  }
172
173
174  size_t matrix::columns(void) const
175  {
176    if (!m_)
177      return 0;
178    return m_->size2;
179  }
180
181
182  gsl_matrix* matrix::create_gsl_matrix_copy(void) const
183  {
184    gsl_matrix* m = gsl_matrix_alloc(rows(),columns());
185    if (!m)
186      throw utility::GSL_error("matrix::create_gsl_matrix_copy failed to allocate memory");
187    if (gsl_matrix_memcpy(m,m_))
188      throw utility::GSL_error("matrix::create_gsl_matrix_copy dimension mis-match");
189    return m;
190  }
191
192
193  int matrix::div_elements(const matrix& b)
194  {
195    return gsl_matrix_div_elements(m_, b.m_);
196  }
197
198
199  bool matrix::equal(const matrix& other, const double d) const
200  {
201    if (columns()!=other.columns() || rows()!=other.rows())
202      return false;
203    for (size_t i=0; i<rows(); i++)
204      for (size_t j=0; j<columns(); j++)
205        // The two last condition checks are needed for NaN detection
206        if (fabs( (*this)(i,j)-other(i,j) ) > d ||
207            (*this)(i,j)!=(*this)(i,j) || other(i,j)!=other(i,j))
208          return false;
209    return true;
210  }
211
212
213  const gsl_matrix* matrix::gsl_matrix_p(void) const
214  {
215    return m_;
216  }
217
218
219  gsl_matrix* matrix::gsl_matrix_p(void)
220  {
221    return m_;
222  }
223
224
225  bool matrix::isnull(void) const
226  {
227    return gsl_matrix_isnull(m_);
228  }
229
230
231  bool matrix::isview(void) const
232  {
233    return view_;
234  }
235
236
237  double matrix::max(void) const
238  {
239    return gsl_matrix_max(m_);
240  }
241
242
243  double matrix::min(void) const
244  {
245    return gsl_matrix_min(m_);
246  }
247
248
249  void matrix::minmax_index(std::pair<size_t,size_t>& min,
250                            std::pair<size_t,size_t>& max) const
251  {
252    gsl_matrix_minmax_index(m_, &min.first, &min.second, &max.first,
253                            &max.second);
254  }
255
256
257  bool matrix::nan(matrix &m) const
258  {
259    m=matrix(rows(),columns(),1.0);
260    bool nan=false;
261    for (size_t i=0; i<rows(); i++)
262      for (size_t j=0; j<columns(); j++) 
263        if (std::isnan(operator()(i,j))) {
264          m(i,j)=0;
265          nan=true;
266        }
267    return nan;
268  }
269
270
271  int matrix::mul_elements(const matrix& b)
272  {
273    return gsl_matrix_mul_elements(m_, b.m_);
274  }
275
276
277  size_t matrix::rows(void) const
278  {
279    if (!m_)
280      return 0;
281    return m_->size1;
282  }
283
284
285  void matrix::scale(const double d)
286  {
287    gsl_matrix_scale(m_, d);
288  }
289
290
291  void matrix::set(const matrix& mat)
292  {
293    if (gsl_matrix_memcpy(m_, mat.m_))
294      throw utility::GSL_error("matrix::create_gsl_matrix_copy dimension mis-match");
295  }
296
297
298  void matrix::set_all(const double value)
299  {
300    gsl_matrix_set_all(m_, value);
301  }
302
303
304  void matrix::set_column(const size_t column, const vector& vec)
305  {
306    int status=gsl_matrix_set_col(m_, column, vec.gsl_vector_p());
307    if (status)
308      throw utility::GSL_error(std::string("matrix::set_column",status));
309  }
310
311
312  void matrix::set_row(const size_t row, const vector& vec)
313  {
314    int status=gsl_matrix_set_row(m_, row, vec.gsl_vector_p());
315    if (status)
316      throw utility::GSL_error(std::string("matrix::set_row",status));
317  }
318
319
320  void matrix::sub(const matrix& b)
321  {
322    int status=gsl_matrix_sub(m_, b.m_);
323    if (status)
324      throw utility::GSL_error(std::string("matrix::sub",status));
325  }
326
327
328  int matrix::swap(matrix& other)
329  {
330    return gsl_matrix_swap(m_, other.m_);
331  }
332
333
334  int matrix::swap_columns(const size_t i, const size_t j)
335  {
336    return gsl_matrix_swap_columns(m_, i, j);
337  }
338
339
340  int matrix::swap_rowcol(const size_t i, const size_t j)
341  {
342    return gsl_matrix_swap_rowcol(m_, i, j);
343  }
344
345
346  int matrix::swap_rows(const size_t i, const size_t j)
347  {
348    return gsl_matrix_swap_rows(m_, i, j);
349  }
350
351
352  void matrix::transpose(void)
353  {
354    if (columns()==rows())
355      gsl_matrix_transpose(m_); // this never fails
356    else {
357      gsl_matrix* transposed = gsl_matrix_alloc(columns(),rows());
358      if (!transposed)
359        throw utility::GSL_error("matrix::transpose failed to allocate memory");
360      // next line never fails if allocation above succeeded.
361      gsl_matrix_transpose_memcpy(transposed,m_);
362      gsl_matrix_free(m_);
363      m_=transposed;
364    }
365  }
366
367
368  double& matrix::operator()(size_t row, size_t column)
369  {
370    return (*gsl_matrix_ptr(m_, row, column));
371  }
372
373
374  const double& matrix::operator()(size_t row, size_t column) const
375  {
376    return (*gsl_matrix_const_ptr(m_, row, column));
377  }
378
379
380  bool matrix::operator==(const matrix& other) const
381  {
382    return equal(other);
383  }
384
385
386  bool matrix::operator!=(const matrix& other) const
387  {
388    return !equal(other);
389  }
390
391
392  const matrix& matrix::operator=( const matrix& other )
393  {
394    if ( this != &other ) {
395      if (view_)
396        delete view_;
397      else if (m_)
398        gsl_matrix_free(m_);
399      m_ = other.create_gsl_matrix_copy();
400    }
401    return *this;
402  }
403
404
405  const matrix& matrix::operator+=(const matrix& m)
406  {
407    add(m);
408    return *this;
409  }
410
411
412  const matrix& matrix::operator+=(const double d)
413  {
414    add(d);
415    return *this;
416  }
417
418
419  const matrix& matrix::operator-=(const matrix& m)
420  {
421    sub(m);
422    return *this;
423  }
424
425
426  const matrix& matrix::operator*=(const matrix& other)
427  {
428    gsl_matrix* result = gsl_matrix_alloc(rows(),other.columns());
429    if (!result)
430      throw utility::GSL_error("matrix::operator*= failed to allocate memory");
431    gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, m_, other.m_, 0.0, result);
432    gsl_matrix_free(m_);
433    m_=result;
434    return *this;
435  }
436
437
438  const matrix& matrix::operator*=(const double d)
439  {
440    scale(d);
441    return *this;
442  }
443
444
445  std::ostream& operator<<(std::ostream& s, const matrix& m)
446  {
447    s.setf(std::ios::dec);
448    s.precision(12);
449    for(size_t i=0, j=0; i<m.rows(); i++)
450      for (j=0; j<m.columns(); j++) {
451        s << m(i,j);
452        if (j<m.columns()-1)
453          s << s.fill();
454        else if (i<m.rows()-1)
455          s << "\n";
456      }
457    return s;
458  }
459
460
461  vector operator*(const matrix& m, const vector& v)
462  {
463    utility::vector res(m.rows());
464    for (size_t i=0; i<res.size(); ++i)
465      res(i) = vector(m,i) * v;
466    return res;
467  }
468
469
470  vector operator*(const vector& v, const matrix& m)
471  {
472    utility::vector res(m.columns());
473    for (size_t i=0; i<res.size(); ++i)
474      res(i) = v * vector(m,i,false);
475    return res;
476  }
477
478}}} // of namespace utility, yat and thep
Note: See TracBrowser for help on using the repository browser.