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

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

Fixes #171 and addresses #2. A few more GSL_error exceptions. Removed Jari comments.

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