source: branches/better_matrix_class/lib/gslapi/matrix.h @ 413

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

All check passes, but more matrix checks needed.

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