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

Last change on this file since 755 was 755, 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.8 KB
Line 
1// $Id: matrix.cc 755 2007-02-18 00:01:39Z 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  void matrix::div_elements(const matrix& b)
194  {
195    int status=gsl_matrix_div_elements(m_, b.m_);
196    if (status)
197      throw utility::GSL_error(std::string("matrix::div_elements",status));
198  }
199
200
201  bool matrix::equal(const matrix& other, const double d) const
202  {
203    if (columns()!=other.columns() || rows()!=other.rows())
204      return false;
205    for (size_t i=0; i<rows(); i++)
206      for (size_t j=0; j<columns(); j++)
207        // The two last condition checks are needed for NaN detection
208        if (fabs( (*this)(i,j)-other(i,j) ) > d ||
209            (*this)(i,j)!=(*this)(i,j) || other(i,j)!=other(i,j))
210          return false;
211    return true;
212  }
213
214
215  const gsl_matrix* matrix::gsl_matrix_p(void) const
216  {
217    return m_;
218  }
219
220
221  gsl_matrix* matrix::gsl_matrix_p(void)
222  {
223    return m_;
224  }
225
226
227  bool matrix::isnull(void) const
228  {
229    return gsl_matrix_isnull(m_);
230  }
231
232
233  bool matrix::isview(void) const
234  {
235    return view_;
236  }
237
238
239  double matrix::max(void) const
240  {
241    return gsl_matrix_max(m_);
242  }
243
244
245  double matrix::min(void) const
246  {
247    return gsl_matrix_min(m_);
248  }
249
250
251  void matrix::minmax_index(std::pair<size_t,size_t>& min,
252                            std::pair<size_t,size_t>& max) const
253  {
254    gsl_matrix_minmax_index(m_, &min.first, &min.second, &max.first,
255                            &max.second);
256  }
257
258
259  bool matrix::nan(matrix &m) const
260  {
261    m=matrix(rows(),columns(),1.0);
262    bool nan=false;
263    for (size_t i=0; i<rows(); i++)
264      for (size_t j=0; j<columns(); j++) 
265        if (std::isnan(operator()(i,j))) {
266          m(i,j)=0;
267          nan=true;
268        }
269    return nan;
270  }
271
272
273  void matrix::mul_elements(const matrix& b)
274  {
275    int status=gsl_matrix_mul_elements(m_, b.m_);
276    if (status)
277      throw utility::GSL_error(std::string("matrix::mul_elements",status));
278   
279  }
280
281
282  size_t matrix::rows(void) const
283  {
284    if (!m_)
285      return 0;
286    return m_->size1;
287  }
288
289
290  void matrix::scale(const double d)
291  {
292    gsl_matrix_scale(m_, d);
293  }
294
295
296  void matrix::set(const matrix& mat)
297  {
298    if (gsl_matrix_memcpy(m_, mat.m_))
299      throw utility::GSL_error("matrix::create_gsl_matrix_copy dimension mis-match");
300  }
301
302
303  void matrix::set_all(const double value)
304  {
305    gsl_matrix_set_all(m_, value);
306  }
307
308
309  void matrix::set_column(const size_t column, const vector& vec)
310  {
311    int status=gsl_matrix_set_col(m_, column, vec.gsl_vector_p());
312    if (status)
313      throw utility::GSL_error(std::string("matrix::set_column",status));
314  }
315
316
317  void matrix::set_row(const size_t row, const vector& vec)
318  {
319    int status=gsl_matrix_set_row(m_, row, vec.gsl_vector_p());
320    if (status)
321      throw utility::GSL_error(std::string("matrix::set_row",status));
322  }
323
324
325  void matrix::sub(const matrix& b)
326  {
327    int status=gsl_matrix_sub(m_, b.m_);
328    if (status)
329      throw utility::GSL_error(std::string("matrix::sub",status));
330  }
331
332
333  void matrix::swap(matrix& other)
334  {
335    int status=gsl_matrix_swap(m_, other.m_);
336    if (status)
337      throw utility::GSL_error(std::string("matrix::swap",status));
338  }
339
340
341  void matrix::swap_columns(const size_t i, const size_t j)
342  {
343    int status=gsl_matrix_swap_columns(m_, i, j);
344    if (status)
345      throw utility::GSL_error(std::string("matrix::swap_columns",status));
346  }
347
348
349  void matrix::swap_rowcol(const size_t i, const size_t j)
350  {
351    int status=gsl_matrix_swap_rowcol(m_, i, j);
352    if (status)
353      throw utility::GSL_error(std::string("matrix::swap_rowcol",status));
354  }
355
356
357  void matrix::swap_rows(const size_t i, const size_t j)
358  {
359    int status=gsl_matrix_swap_rows(m_, i, j);
360    if (status)
361      throw utility::GSL_error(std::string("matrix::swap_rows",status));
362  }
363
364
365  void matrix::transpose(void)
366  {
367    if (columns()==rows())
368      gsl_matrix_transpose(m_); // this never fails
369    else {
370      gsl_matrix* transposed = gsl_matrix_alloc(columns(),rows());
371      if (!transposed)
372        throw utility::GSL_error("matrix::transpose failed to allocate memory");
373      // next line never fails if allocation above succeeded.
374      gsl_matrix_transpose_memcpy(transposed,m_);
375      gsl_matrix_free(m_);
376      m_=transposed;
377    }
378  }
379
380
381  double& matrix::operator()(size_t row, size_t column)
382  {
383    double* d=gsl_matrix_ptr(m_, row, column);
384    if (!d)
385      throw utility::GSL_error("matrix::operator()",GSL_EINVAL);
386    return *d;
387  }
388
389
390  const double& matrix::operator()(size_t row, size_t column) const
391  {
392    const double* d=gsl_matrix_const_ptr(m_, row, column);
393    if (!d)
394      throw utility::GSL_error("matrix::operator()",GSL_EINVAL);
395    return *d;
396  }
397
398
399  bool matrix::operator==(const matrix& other) const
400  {
401    return equal(other);
402  }
403
404
405  bool matrix::operator!=(const matrix& other) const
406  {
407    return !equal(other);
408  }
409
410
411  const matrix& matrix::operator=( const matrix& other )
412  {
413    if ( this != &other ) {
414      if (view_)
415        delete view_;
416      else if (m_)
417        gsl_matrix_free(m_);
418      m_ = other.create_gsl_matrix_copy();
419    }
420    return *this;
421  }
422
423
424  const matrix& matrix::operator+=(const matrix& m)
425  {
426    add(m);
427    return *this;
428  }
429
430
431  const matrix& matrix::operator+=(const double d)
432  {
433    add(d);
434    return *this;
435  }
436
437
438  const matrix& matrix::operator-=(const matrix& m)
439  {
440    sub(m);
441    return *this;
442  }
443
444
445  const matrix& matrix::operator*=(const matrix& other)
446  {
447    gsl_matrix* result = gsl_matrix_alloc(rows(),other.columns());
448    if (!result)
449      throw utility::GSL_error("matrix::operator*= failed to allocate memory");
450    gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, m_, other.m_, 0.0, result);
451    gsl_matrix_free(m_);
452    m_=result;
453    return *this;
454  }
455
456
457  const matrix& matrix::operator*=(const double d)
458  {
459    scale(d);
460    return *this;
461  }
462
463
464  std::ostream& operator<<(std::ostream& s, const matrix& m)
465  {
466    s.setf(std::ios::dec);
467    s.precision(12);
468    for(size_t i=0, j=0; i<m.rows(); i++)
469      for (j=0; j<m.columns(); j++) {
470        s << m(i,j);
471        if (j<m.columns()-1)
472          s << s.fill();
473        else if (i<m.rows()-1)
474          s << "\n";
475      }
476    return s;
477  }
478
479
480  vector operator*(const matrix& m, const vector& v)
481  {
482    utility::vector res(m.rows());
483    for (size_t i=0; i<res.size(); ++i)
484      res(i) = vector(m,i) * v;
485    return res;
486  }
487
488
489  vector operator*(const vector& v, const matrix& m)
490  {
491    utility::vector res(m.columns());
492    for (size_t i=0; i<res.size(); ++i)
493      res(i) = v * vector(m,i,false);
494    return res;
495  }
496
497}}} // of namespace utility, yat and thep
Note: See TracBrowser for help on using the repository browser.