SourceForge Logo Tiny Vector Matrix library using Expression Templates Sourceforge Project Page

include/tvmet/Matrix.h

Go to the documentation of this file.
00001 /*
00002  * Tiny Vector Matrix Library
00003  * Dense Vector Matrix Libary of Tiny size using Expression Templates
00004  *
00005  * Copyright (C) 2001 - 2007 Olaf Petzold <opetzold@users.sourceforge.net>
00006  *
00007  * This library is free software; you can redistribute it and/or
00008  * modify it under the terms of the GNU lesser General Public
00009  * License as published by the Free Software Foundation; either
00010  * version 2.1 of the License, or (at your option) any later version.
00011  *
00012  * This library is distributed in the hope that it will be useful,
00013  * but WITHOUT ANY WARRANTY; without even the implied warranty of
00014  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00015  * lesser General Public License for more details.
00016  *
00017  * You should have received a copy of the GNU lesser General Public
00018  * License along with this library; if not, write to the Free Software
00019  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00020  *
00021  * $Id: Matrix.h,v 1.58 2007-06-23 15:58:58 opetzold Exp $
00022  */
00023 
00024 #ifndef TVMET_MATRIX_H
00025 #define TVMET_MATRIX_H
00026 
00027 #include <iterator>         // reverse_iterator
00028 
00029 #include <tvmet/tvmet.h>
00030 #include <tvmet/TypePromotion.h>
00031 #include <tvmet/CommaInitializer.h>
00032 #include <tvmet/RunTimeError.h>
00033 
00034 #include <tvmet/xpr/Matrix.h>
00035 #include <tvmet/xpr/MatrixRow.h>
00036 #include <tvmet/xpr/MatrixCol.h>
00037 #include <tvmet/xpr/MatrixDiag.h>
00038 
00039 namespace tvmet {
00040 
00041 
00042 /* forwards */
00043 template<class T, std::size_t Rows, std::size_t Cols> class Matrix;
00044 template<class T,
00045    std::size_t RowsBgn, std::size_t RowsEnd,
00046    std::size_t ColsBgn, std::size_t ColsEnd,
00047    std::size_t RowStride, std::size_t ColStride /*=1*/>
00048 class MatrixSliceConstReference; // unused here; for me only
00049 
00050 
00055 template<class T, std::size_t NRows, std::size_t NCols>
00056 class MatrixConstReference
00057   : public TvmetBase < MatrixConstReference<T, NRows, NCols> >
00058 {
00059 public:
00060   typedef T           value_type;
00061   typedef T*            pointer;
00062   typedef const T*          const_pointer;
00063 
00065   enum {
00066     Rows = NRows,     
00067     Cols = NCols,     
00068     Size = Rows * Cols      
00069   };
00070 
00071 public:
00073   enum {
00074     ops       = Rows * Cols
00075   };
00076 
00077 private:
00078   MatrixConstReference();
00079   MatrixConstReference& operator=(const MatrixConstReference&);
00080 
00081 public:
00083   explicit MatrixConstReference(const Matrix<T, Rows, Cols>& rhs)
00084     : m_data(rhs.data())
00085   { }
00086 
00088   explicit MatrixConstReference(const_pointer data)
00089     : m_data(data)
00090   { }
00091 
00092 public: // access operators
00094   value_type operator()(std::size_t i, std::size_t j) const {
00095     TVMET_RT_CONDITION((i < Rows) && (j < Cols), "MatrixConstReference Bounce Violation")
00096     return m_data[i * Cols + j];
00097   }
00098 
00099 public: // debugging Xpr parse tree
00100   void print_xpr(std::ostream& os, std::size_t l=0) const {
00101     os << IndentLevel(l)
00102        << "MatrixConstReference[O=" << ops << "]<"
00103        << "T=" << typeid(value_type).name() << ">,"
00104        << std::endl;
00105   }
00106 
00107 private:
00108   const_pointer _tvmet_restrict       m_data;
00109 };
00110 
00111 
00122 template<class T, std::size_t NRows, std::size_t NCols>
00123 class Matrix
00124 {
00125 public:
00127   typedef T           value_type;
00128 
00130   typedef T&              reference;
00131 
00133   typedef const T&              const_reference;
00134 
00136   typedef T*              iterator;
00137 
00139   typedef const T*              const_iterator;
00140 
00142   typedef std::reverse_iterator<iterator>     reverse_iterator;
00143 
00145   typedef std::reverse_iterator<const_iterator>   const_reverse_iterator;
00146 
00147 public:
00149   enum {
00150     Rows = NRows,     
00151     Cols = NCols,     
00152     Size = Rows * Cols      
00153   };
00154 
00155 public:
00157   enum {
00158     ops_assign = Rows * Cols,
00159     ops        = ops_assign,
00160     use_meta   = ops < TVMET_COMPLEXITY_M_ASSIGN_TRIGGER ? true : false
00161   };
00162 
00163 public: // STL  interface
00165   iterator begin() { return m_data; }
00166 
00168   iterator end() { return m_data + Size; }
00169 
00171   const_iterator begin() const { return m_data; }
00172 
00174   const_iterator end() const { return m_data + Size; }
00175 
00177   reverse_iterator rbegin() { return reverse_iterator( end() ); }
00178 
00180   const_reverse_iterator rbegin() const {
00181     return const_reverse_iterator( end() );
00182   }
00183 
00185   reverse_iterator rend() { return reverse_iterator( begin() ); }
00186 
00188   const_reverse_iterator rend() const {
00189     return const_reverse_iterator( begin() );
00190   }
00191 
00193   static std::size_t size() { return Size; }
00194 
00196   static std::size_t max_size() { return Size; }
00197 
00199   static bool empty() { return false; }
00200 
00201 public:
00203   static std::size_t rows() { return Rows; }
00204 
00206   static std::size_t cols() { return Cols; }
00207 
00208 public:
00210   ~Matrix() {
00211 #if defined(TVMET_DYNAMIC_MEMORY)
00212     delete [] m_data;
00213 #endif
00214   }
00215 
00218   explicit Matrix()
00219 #if defined(TVMET_DYNAMIC_MEMORY)
00220     : m_data( new value_type[Size] )
00221 #endif
00222   { }
00223 
00225   Matrix(const Matrix& rhs)
00226 #if defined(TVMET_DYNAMIC_MEMORY)
00227     : m_data( new value_type[Size] )
00228 #endif
00229   {
00230     *this = XprMatrix<ConstReference, Rows, Cols>(rhs.const_ref());
00231   }
00232 
00237   template<class InputIterator>
00238   explicit Matrix(InputIterator first, InputIterator last)
00239 #if defined(TVMET_DYNAMIC_MEMORY)
00240     : m_data( new value_type[Size] )
00241 #endif
00242   {
00243     TVMET_RT_CONDITION(static_cast<std::size_t>(std::distance(first, last)) <= Size,
00244            "InputIterator doesn't fits in size" )
00245     std::copy(first, last, m_data);
00246   }
00247 
00252   template<class InputIterator>
00253   explicit Matrix(InputIterator first, std::size_t sz)
00254 #if defined(TVMET_DYNAMIC_MEMORY)
00255     : m_data( new value_type[Size] )
00256 #endif
00257   {
00258     TVMET_RT_CONDITION(sz <= Size, "InputIterator doesn't fits in size" )
00259     std::copy(first, first + sz, m_data);
00260   }
00261 
00263   explicit Matrix(value_type rhs)
00264 #if defined(TVMET_DYNAMIC_MEMORY)
00265     : m_data( new value_type[Size] )
00266 #endif
00267   {
00268     typedef XprLiteral<value_type> expr_type;
00269     *this = XprMatrix<expr_type, Rows, Cols>(expr_type(rhs));
00270   }
00271 
00273   template<class E>
00274   explicit Matrix(const XprMatrix<E, Rows, Cols>& e)
00275 #if defined(TVMET_DYNAMIC_MEMORY)
00276     : m_data( new value_type[Size] )
00277 #endif
00278   {
00279     *this = e;
00280   }
00281 
00284   CommaInitializer<Matrix, Size> operator=(value_type rhs) {
00285     return CommaInitializer<Matrix, Size>(*this, rhs);
00286   }
00287 
00288 public: // access operators
00289   value_type* _tvmet_restrict data() { return m_data; }
00290   const value_type* _tvmet_restrict data() const { return m_data; }
00291 
00292 public: // index access operators
00293   value_type& _tvmet_restrict operator()(std::size_t i, std::size_t j) {
00294     // Note: g++-2.95.3 does have problems on typedef reference
00295     TVMET_RT_CONDITION((i < Rows) && (j < Cols), "Matrix Bounce Violation")
00296     return m_data[i * Cols + j];
00297   }
00298 
00299   value_type operator()(std::size_t i, std::size_t j) const {
00300     TVMET_RT_CONDITION((i < Rows) && (j < Cols), "Matrix Bounce Violation")
00301     return m_data[i * Cols + j];
00302   }
00303 
00304 public: // ET interface
00305   typedef MatrixConstReference<T, Rows, Cols>     ConstReference;
00306 
00307   typedef MatrixSliceConstReference<
00308     T,
00309     0, Rows, 0, Cols,
00310     Rows, 1
00311   >             SliceConstReference;
00312 
00314   ConstReference const_ref() const { return ConstReference(*this); }
00315 
00320   ConstReference const_sliceref() const { return SliceConstReference(*this); }
00321 
00323   XprMatrix<ConstReference, Rows, Cols> as_expr() const {
00324     return XprMatrix<ConstReference, Rows, Cols>(this->const_ref());
00325   }
00326 
00327 private:
00329   template<class Dest, class Src, class Assign>
00330   static inline
00331   void do_assign(dispatch<true>, Dest& dest, const Src& src, const Assign& assign_fn) {
00332     meta::Matrix<Rows, Cols, 0, 0>::assign(dest, src, assign_fn);
00333   }
00334 
00336   template<class Dest, class Src, class Assign>
00337   static inline
00338   void do_assign(dispatch<false>, Dest& dest, const Src& src, const Assign& assign_fn) {
00339     loop::Matrix<Rows, Cols>::assign(dest, src, assign_fn);
00340   }
00341 
00342 private:
00345   template<class T2, class Assign>
00346   void assign_to(Matrix<T2, Rows, Cols>& dest, const Assign& assign_fn) const {
00347     do_assign(dispatch<use_meta>(), dest, *this, assign_fn);
00348   }
00349 
00350 public:  // assign operations
00354   template<class T2>
00355   Matrix& operator=(const Matrix<T2, Rows, Cols>& rhs) {
00356     rhs.assign_to(*this, Fcnl_assign<value_type, T2>());
00357     return *this;
00358   }
00359 
00361   template <class E>
00362   Matrix& operator=(const XprMatrix<E, Rows, Cols>& rhs) {
00363     rhs.assign_to(*this, Fcnl_assign<value_type, typename E::value_type>());
00364     return *this;
00365   }
00366 
00367 private:
00368   template<class Obj, std::size_t LEN> friend class CommaInitializer;
00369 
00373   Matrix& assign_value(value_type rhs) {
00374     typedef XprLiteral<value_type>      expr_type;
00375     *this = XprMatrix<expr_type, Rows, Cols>(expr_type(rhs));
00376     return *this;
00377   }
00378 
00379 public: // math operators with scalars
00380   // NOTE: this meaning is clear - element wise ops even if not in ns element_wise
00381   Matrix& operator+=(value_type) TVMET_CXX_ALWAYS_INLINE;
00382   Matrix& operator-=(value_type) TVMET_CXX_ALWAYS_INLINE;
00383   Matrix& operator*=(value_type) TVMET_CXX_ALWAYS_INLINE;
00384   Matrix& operator/=(value_type) TVMET_CXX_ALWAYS_INLINE;
00385 
00386   Matrix& operator%=(std::size_t) TVMET_CXX_ALWAYS_INLINE;
00387   Matrix& operator^=(std::size_t) TVMET_CXX_ALWAYS_INLINE;
00388   Matrix& operator&=(std::size_t) TVMET_CXX_ALWAYS_INLINE;
00389   Matrix& operator|=(std::size_t) TVMET_CXX_ALWAYS_INLINE;
00390   Matrix& operator<<=(std::size_t) TVMET_CXX_ALWAYS_INLINE;
00391   Matrix& operator>>=(std::size_t) TVMET_CXX_ALWAYS_INLINE;
00392 
00393 public: // math operators with matrizes
00394   // NOTE: access using the operators in ns element_wise, since that's what is does
00395   template <class T2> Matrix& M_add_eq(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00396   template <class T2> Matrix& M_sub_eq(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00397   template <class T2> Matrix& M_mul_eq(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00398   template <class T2> Matrix& M_div_eq(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00399   template <class T2> Matrix& M_mod_eq(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00400   template <class T2> Matrix& M_xor_eq(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00401   template <class T2> Matrix& M_and_eq(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00402   template <class T2> Matrix& M_or_eq (const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00403   template <class T2> Matrix& M_shl_eq(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00404   template <class T2> Matrix& M_shr_eq(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00405 
00406 public: // math operators with expressions
00407   // NOTE: access using the operators in ns element_wise, since that's what is does
00408   template <class E> Matrix& M_add_eq(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00409   template <class E> Matrix& M_sub_eq(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00410   template <class E> Matrix& M_mul_eq(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00411   template <class E> Matrix& M_div_eq(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00412   template <class E> Matrix& M_mod_eq(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00413   template <class E> Matrix& M_xor_eq(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00414   template <class E> Matrix& M_and_eq(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00415   template <class E> Matrix& M_or_eq (const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00416   template <class E> Matrix& M_shl_eq(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00417   template <class E> Matrix& M_shr_eq(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00418 
00419 public: // aliased math operators with expressions
00420   template <class T2> Matrix& alias_assign(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00421   template <class T2> Matrix& alias_add_eq(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00422   template <class T2> Matrix& alias_sub_eq(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00423   template <class T2> Matrix& alias_mul_eq(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00424   template <class T2> Matrix& alias_div_eq(const Matrix<T2, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00425 
00426   template <class E> Matrix& alias_assign(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00427   template <class E> Matrix& alias_add_eq(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00428   template <class E> Matrix& alias_sub_eq(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00429   template <class E> Matrix& alias_mul_eq(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00430   template <class E> Matrix& alias_div_eq(const XprMatrix<E, Rows, Cols>&) TVMET_CXX_ALWAYS_INLINE;
00431 
00432 public: // io
00434   struct Info : public TvmetBase<Info> {
00435     std::ostream& print_xpr(std::ostream& os) const {
00436       os << "Matrix<T=" << typeid(value_type).name()
00437    << ", R=" << Rows << ", C=" << Cols << ">";
00438       return os;
00439     }
00440   };
00441 
00443   static Info info() { return Info(); }
00444 
00446   std::ostream& print_xpr(std::ostream& os, std::size_t l=0) const;
00447 
00449   std::ostream& print_on(std::ostream& os) const;
00450 
00451 private:
00453 #if defined(TVMET_DYNAMIC_MEMORY)
00454   value_type*           m_data;
00455 #else
00456   value_type            m_data[Size];
00457 #endif
00458 };
00459 
00460 
00461 } // namespace tvmet
00462 
00463 #include <tvmet/MatrixImpl.h>
00464 #include <tvmet/MatrixFunctions.h>
00465 #include <tvmet/MatrixBinaryFunctions.h>
00466 #include <tvmet/MatrixUnaryFunctions.h>
00467 #include <tvmet/MatrixOperators.h>
00468 #include <tvmet/MatrixEval.h>
00469 #include <tvmet/AliasProxy.h>
00470 
00471 #endif // TVMET_MATRIX_H
00472 
00473 // Local Variables:
00474 // mode:C++
00475 // tab-width:8
00476 // End:

Author: