source: trunk/yat/utility/matrix.h @ 703

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

Addresses #65 and #170.

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 13.6 KB
Line 
1#ifndef _theplu_yat_utility_matrix_
2#define _theplu_yat_utility_matrix_
3
4// $Id: matrix.h 703 2006-12-18 00:47:44Z jari $
5
6/*
7  Copyright (C) 2003 Daniel Dalevi, Peter Johansson
8  Copyright (C) 2004 Jari Häkkinen, Peter Johansson
9  Copyright (C) 2005 Jari Häkkinen, Peter Johansson, Markus Ringnér
10  Copyright (C) 2006 Jari Häkkinen, Peter Johansson
11
12  This file is part of the yat library, http://lev.thep.lu.se/trac/yat
13
14  The yat library is free software; you can redistribute it and/or
15  modify it under the terms of the GNU General Public License as
16  published by the Free Software Foundation; either version 2 of the
17  License, or (at your option) any later version.
18
19  The yat library is distributed in the hope that it will be useful,
20  but WITHOUT ANY WARRANTY; without even the implied warranty of
21  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
22  General Public License for more details.
23
24  You should have received a copy of the GNU General Public License
25  along with this program; if not, write to the Free Software
26  Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA
27  02111-1307, USA.
28*/
29
30#include "vector.h"
31#include "Exception.h"
32
33#include <gsl/gsl_matrix.h>
34#include <iostream>
35
36namespace theplu {
37namespace yat {
38namespace utility {
39
40  class vector;
41
42  ///
43  /// This is the yat interface to GSL matrix. 'double' is the
44  /// only type supported, maybe we should add a 'complex' type as
45  /// well in the future.
46  ///
47  /// \par[File streams] Reading and writing vectors to file streams
48  /// are of course supported. These are implemented without using GSL
49  /// functionality, and thus binary read and write to streams are not
50  /// supported.
51  ///
52  /// \par[Matrix views] GSL matrix views are supported and these are
53  /// disguised as ordinary utility::matrix'es. A support function is
54  /// added, utility::matrix::isview(), that can be used to check if a
55  /// matrix object is a view. Note that view matricess do not own the
56  /// undelying data, and a view is not valid if the matrix owning the
57  /// data is deallocated.
58  ///
59  /// @note All GSL matrix related functions are not implement but
60  /// most functionality defined for GSL matrices can be achieved with
61  /// this interface class. Most notable GSL functionality not
62  /// supported are no binary file support and views on arrays,
63  /// utility::vectors, gsl_vectors, diagonals, subdiagonals, and
64  /// superdiagonals. If there is an interest from the user community,
65  /// the omitted functionality could be included.
66  ///
67  /// @todo Maybe it would be smart to create temporary objects need
68  /// for BLAS in constructors. As it is now temporary objects are
69  /// called before BLAS functionality iss used, cf. const matrix&
70  /// matrix::operator*=(const matrix& other)
71  ///
72  class matrix
73  {
74  public:
75
76    /**
77       @brief The default constructor.
78   
79       This contructor does not initialize underlying (essential)
80       structures.
81    */
82    matrix(void);
83
84    ///
85    /// @brief Constructor allocating memory space for \a r times \a c
86    /// elements, and sets all elements to \a init_value.
87    ///
88    matrix(const size_t& r, const size_t& c, double init_value=0);
89
90    ///
91    /// @brief The copy constructor.
92    ///
93    /// @note If the object to be copied is a matrix view, the values
94    /// of the view will be copied, i.e. the view is not copied.
95    ///
96    matrix(const matrix&);
97
98    ///
99    /// The matrix view constructor.
100    ///
101    /// Create a view of matrix \a m, with starting row \a offset_row,
102    /// starting column \a offset_column, row size \a n_row, and
103    /// column size \a n_column.
104    ///
105    /// A matrix view can be used as any matrix with the difference
106    /// that changes made to the view will also change the object that
107    /// is viewed. Also, using the copy constructor will create a new
108    /// matrix object that is a copy of whatever is viewed. If a copy
109    /// of the view is needed then you should use this constructor to
110    /// obtain a copy.
111    ///
112    /// @note If the object viewed by the view goes out of scope or is
113    /// deleted, the view becomes invalid and the result of further
114    /// use is undefined.
115    ///
116    matrix(matrix& m, size_t offset_row, size_t offset_column, size_t n_row,
117           size_t n_column);
118
119    ///
120    /// The istream constructor.
121    ///
122    /// Matrix rows are sepearated with the new line character, and
123    /// column elements in a row must be separated either with white space
124    /// characters except the new line (\\n) character (default), or
125    /// by the delimiter \a sep.
126    ///. Rows, and
127    /// columns, are read consecutively and only rectangular matrices
128    /// are supported. Empty lines are ignored. End of input is at end
129    /// of file marker.
130    ///
131    explicit matrix(std::istream &, char sep='\0') 
132      throw(utility::IO_error, std::exception);
133
134    ///
135    /// @brief The destructor.
136    ///
137    ~matrix(void);
138
139    ///
140    /// Elementwise addition of the elements of matrix \a b to the
141    /// elements of the calling matrix ,\f$ a_{ij} = a_{ij} + b_{ij} \;
142    /// \forall i,j \f$. The result is stored into the calling matrix.
143    ///
144    /// @return Whatever GSL returns.
145    ///
146    inline int add(const matrix& b) { return gsl_matrix_add(m_,b.m_); }
147
148    ///
149    /// Add the scalar value \a d to the elements of the calling
150    /// matrix, \f$ a_{ij} = a_{ij} + d \; \forall i,j \f$. The result
151    /// is stored into the calling matrix.
152    ///
153    /// @return Whatever GSL returns.
154    ///
155    inline int
156    add_constant(const double d) { return gsl_matrix_add_constant(m_,d); }
157
158    ///
159    /// @return The number of columns in the matrix.
160    ///
161    inline size_t columns(void) const { return m_->size2; }
162
163    ///
164    /// Elementwise division of the elemnts of the calling matrix by
165    /// the elements of matrix \a b, \f$ a_{ij} = a_{ij} / b_{ij} \;
166    /// \forall i,j \f$. The result is stored into the calling matrix.
167    ///
168    /// @return Whatever GSL returns.
169    ///
170    inline int
171    div_elements(const matrix& b) { return gsl_matrix_div_elements(m_,b.m_); }
172
173    ///
174    /// Check whether matrices are equal within a user defined
175    /// precision, set by \a precision.
176    ///
177    /// @return True if each element deviates less or equal than \a
178    /// d. If any matrix contain a NaN, false is always returned.
179    ///
180    bool equal(const matrix&, const double precision=0) const;
181
182    ///
183    /// @return A const pointer to the internal GSL matrix.
184    ///
185    inline const gsl_matrix* gsl_matrix_p(void) const { return m_; }
186
187    ///
188    /// @return A pointer to the internal GSL matrix.
189    ///
190    // Jari, is this needed?
191    inline gsl_matrix* gsl_matrix_p(void) { return m_; };
192
193    ///
194    /// @return True if all elements in the matrix is zero, false
195    /// othwerwise;
196    ///
197    inline bool isnull(void) const { return gsl_matrix_isnull(m_); }
198
199    ///
200    /// @brief Check if the matrix object is a view (sub-matrix) to
201    /// another matrix.
202    ///
203    /// @return True if the object is a view, false othwerwise.
204    ///
205    inline bool isview(void) const { return view_; }
206
207    ///
208    /// @return The maximum value of the matrix.
209    ///
210    inline double max(void) const { return gsl_matrix_max(m_); }
211
212    ///
213    /// @return The index to the element with the maximum value of the
214    /// matrix. The lowest index has precedence (searching in
215    /// row-major order).
216    ///
217    std::pair<size_t,size_t> max_index(void) const;
218
219    ///
220    /// @return The minimum value of the matrix.
221    ///
222    inline double min(void) const { return gsl_matrix_min(m_); }
223
224    ///
225    /// @return The index to the element with the minimum value of the
226    /// matrix. The lowest index has precedence (searching in
227    /// row-major order).
228    ///
229    std::pair<size_t,size_t> min_index(void) const;
230
231    ///
232    /// @return The indecies to the element with the minimum and
233    /// maximum values of the matrix, respectively. The lowest index
234    /// has precedence (searching in row-major order).
235    ///
236    inline void minmax_index(std::pair<size_t,size_t>& min,
237                             std::pair<size_t,size_t>& max) const
238      { gsl_matrix_minmax_index(m_,&min.first,&min.second,
239                                &max.first,&max.second); }
240
241    ///
242    /// @return The minimum and maximum values of the matrix, as the
243    /// \a first and \a second member of the returned \a pair,
244    /// respectively.
245    ///
246    std::pair<double,double> minmax(void) const;
247
248    ///
249    /// The argument is changed into a matrix containing 1's and 0's
250    /// only. 1 means element in matrix is valid and a zero means
251    /// element is a NaN.
252    ///
253    /// @return true if matrix contains at least one NaN.
254    ///
255    bool nan(matrix&) const;
256
257    ///
258    /// Multiply the elements of matrix \a b with the elements of the
259    /// calling matrix ,\f$ a_{ij} = a_{ij} * b_{ij} \; \forall
260    /// i,j \f$. The result is stored into the calling matrix.
261    ///
262    /// @return Whatever GSL returns.
263    ///
264    inline int
265    mul_elements(const matrix& b) { return gsl_matrix_mul_elements(m_,b.m_); }
266
267    ///
268    /// @return The number of rows in the matrix.
269    ///
270    inline size_t rows(void) const { return m_->size1; }
271
272    ///
273    /// Multiply the elements of the calling matrix with a scalar \a
274    /// d, \f$ a_{ij} = d * a_{ij} \; \forall i,j \f$. The result is
275    /// stored into the calling matrix.
276    ///
277    /// @return Whatever GSL returns.
278    ///
279    inline int scale(const double d) { return gsl_matrix_scale(m_,d); }
280
281    ///
282    /// Set element values to values in \a mat. This function is
283    /// needed for assignment of viewed elements.
284    ///
285    /// @return Whatever GSL returns.
286    ///
287    /// @note No check on size is done.
288    ///
289    /// @see const matrix& operator=(const matrix&)
290    ///
291    inline int set(const matrix& mat) { return gsl_matrix_memcpy(m_,mat.m_); }
292
293    ///
294    /// Set all elements to \a value.
295    ///
296    inline void set_all(const double value) { gsl_matrix_set_all(m_,value); }
297
298    ///
299    /// Set \a column values to values in \a vec.
300    ///
301    /// @return Whatever GSL returns.
302    ///
303    /// @note No check on size is done.
304    ///
305    inline int set_column(const size_t column, const vector& vec)
306      { return gsl_matrix_set_col(m_, column, vec.gsl_vector_p()); }
307
308    ///
309    /// @brief Set \a row values to values in \a vec.
310    ///
311    /// @return Whatever GSL returns.
312    ///
313    /// @note No check on size is done.
314    ///
315    inline int set_row(const size_t row, const vector& vec)
316      { return gsl_matrix_set_row(m_, row, vec.gsl_vector_p()); }
317
318    ///
319    /// Subtract the elements of matrix \a b from the elements of the
320    /// calling matrix ,\f$ a_{ij} = a_{ij} - b_{ij} \; \forall
321    /// i,j \f$. The result is stored into the calling matrix.
322    ///
323    /// @return Whatever GSL returns.
324    ///
325    inline int sub(const matrix& b) { return gsl_matrix_sub(m_,b.m_); }
326
327    ///
328    /// Exchange the elements of the this and \a other by copying. The
329    /// two matrices must have the same size.
330    ///
331    /// @return Whatever GSL returns.
332    ///
333    inline int swap(matrix& other) { return gsl_matrix_swap(m_,other.m_); }
334
335    ///
336    /// @brief Swap columns \a i and \a j.
337    ///
338    inline int swap_columns(const size_t i,const size_t j)
339      { return gsl_matrix_swap_columns(m_,i,j); }
340
341    ///
342    /// @brief Swap row \a i and column \a j.
343    ///
344    inline int swap_rowcol(const size_t i,const size_t j)
345      { return gsl_matrix_swap_rowcol(m_,i,j); }
346
347    ///
348    /// @brief Swap rows \a i and \a j.
349    ///
350    inline int swap_rows(const size_t i, const size_t j)
351      { return gsl_matrix_swap_rows(m_,i,j); }
352
353    ///
354    /// @brief Transpose the matrix.
355    ///
356    void transpose(void);
357
358    ///
359    /// @return Reference to the element position (\a row, \a column).
360    ///
361    inline double& operator()(size_t row,size_t column)
362    { return (*gsl_matrix_ptr(m_,row,column)); }
363
364    ///
365    /// @return Const reference to the element position (\a row, \a
366    /// column).
367    ///
368    inline const double& operator()(size_t row,size_t column) const
369    { return (*gsl_matrix_const_ptr(m_,row,column)); }
370
371    ///
372    /// @brief Matrix-vector multiplication.
373    ///
374    /// @return The resulting vector.
375    ///
376    // Jari, where should this go?
377    const vector operator*(const vector&) const;
378
379    ///
380    /// @brief Comparison operator.
381    ///
382    /// @return True if all elements are equal otherwise False.
383    ///
384    /// @see equal
385    ///
386    inline bool operator==(const matrix& other) const { return equal(other); }
387
388    ///
389    /// @brief Comparison operator.
390    ///
391    /// @return False if all elements are equal otherwise True.
392    ///
393    /// @see equal
394    ///
395    inline bool operator!=(const matrix& other) const { return !equal(other); }
396
397    ///
398    /// @brief The assignment operator.
399    ///
400    /// There is no requirements on dimensions, i.e. the matrix is
401    /// remapped in memory if necessary. This implies that in general
402    /// views cannot be assigned using this operator. Views will be
403    /// mutated into normal matrices. The only exception to this
404    /// behaviour on views is when self-assignemnt is done, since
405    /// self-assignment is ignored.
406    ///
407    /// @return A const reference to the resulting matrix.
408    ///
409    /// @see int set(const matrix&)
410    ///
411    const matrix& operator=(const matrix& other);
412
413    ///
414    /// @brief Add and assign operator.
415    ///
416    const matrix& operator+=(const matrix&);
417
418    ///
419    /// @brief Add and assign operator.
420    ///
421    const matrix& operator+=(const double);
422
423    ///
424    /// @brief Subtract and assign operator.
425    ///
426    const matrix& operator-=(const matrix&);
427
428    ///
429    /// @brief Multiply and assigment operator.
430    ///
431    /// @return Const reference to the resulting matrix.
432    ///
433    const matrix& operator*=(const matrix&);
434
435    ///
436    /// @brief Multiply and assignment operator
437    ///
438    const matrix& operator*=(const double);
439
440  private:
441
442    ///
443    /// Create a new copy of the internal GSL matrix.
444    ///
445    /// Necessary memory for the new GSL matrix is allocated and the
446    /// caller is responsible for freeing the allocated memory.
447    ///
448    /// @return A pointer to a copy of the internal GSL matrix.
449    ///
450    gsl_matrix* create_gsl_matrix_copy(void) const;
451
452    gsl_matrix* m_;
453    gsl_matrix_view* view_;
454  };
455
456  ///
457  /// The output operator for the matrix class.
458  ///
459  std::ostream& operator<< (std::ostream& s, const matrix&);
460
461
462}}} // of namespace utility, yat, and theplu
463
464#endif
Note: See TracBrowser for help on using the repository browser.