plugins/base1/se.lu.thep.wenni/trunk/lib/c++_tools/gslapi/matrix.h

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