eigen3_interface.hh
Go to the documentation of this file.
00001 //=====================================================
00002 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
00003 //=====================================================
00004 //
00005 // This program is free software; you can redistribute it and/or
00006 // modify it under the terms of the GNU General Public License
00007 // as published by the Free Software Foundation; either version 2
00008 // of the License, or (at your option) any later version.
00009 //
00010 // This program is distributed in the hope that it will be useful,
00011 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00012 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00013 // GNU General Public License for more details.
00014 // You should have received a copy of the GNU General Public License
00015 // along with this program; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA  02111-1307, USA.
00017 //
00018 #ifndef EIGEN3_INTERFACE_HH
00019 #define EIGEN3_INTERFACE_HH
00020 
00021 #include <Eigen/Eigen>
00022 #include <vector>
00023 #include "btl.hh"
00024 
00025 using namespace Eigen;
00026 
00027 template<class real, int SIZE=Dynamic>
00028 class eigen3_interface
00029 {
00030 
00031 public :
00032 
00033   enum {IsFixedSize = (SIZE!=Dynamic)};
00034 
00035   typedef real real_type;
00036 
00037   typedef std::vector<real> stl_vector;
00038   typedef std::vector<stl_vector> stl_matrix;
00039 
00040   typedef Eigen::Matrix<real,SIZE,SIZE> gene_matrix;
00041   typedef Eigen::Matrix<real,SIZE,1> gene_vector;
00042 
00043   static inline std::string name( void )
00044   {
00045     return EIGEN_MAKESTRING(BTL_PREFIX);
00046   }
00047 
00048   static void free_matrix(gene_matrix & A, int N) {}
00049 
00050   static void free_vector(gene_vector & B) {}
00051 
00052   static BTL_DONT_INLINE void matrix_from_stl(gene_matrix & A, stl_matrix & A_stl){
00053     A.resize(A_stl[0].size(), A_stl.size());
00054 
00055     for (int j=0; j<A_stl.size() ; j++){
00056       for (int i=0; i<A_stl[j].size() ; i++){
00057         A.coeffRef(i,j) = A_stl[j][i];
00058       }
00059     }
00060   }
00061 
00062   static BTL_DONT_INLINE  void vector_from_stl(gene_vector & B, stl_vector & B_stl){
00063     B.resize(B_stl.size(),1);
00064 
00065     for (int i=0; i<B_stl.size() ; i++){
00066       B.coeffRef(i) = B_stl[i];
00067     }
00068   }
00069 
00070   static BTL_DONT_INLINE  void vector_to_stl(gene_vector & B, stl_vector & B_stl){
00071     for (int i=0; i<B_stl.size() ; i++){
00072       B_stl[i] = B.coeff(i);
00073     }
00074   }
00075 
00076   static BTL_DONT_INLINE  void matrix_to_stl(gene_matrix & A, stl_matrix & A_stl){
00077     int N=A_stl.size();
00078 
00079     for (int j=0;j<N;j++){
00080       A_stl[j].resize(N);
00081       for (int i=0;i<N;i++){
00082         A_stl[j][i] = A.coeff(i,j);
00083       }
00084     }
00085   }
00086 
00087   static inline void matrix_matrix_product(const gene_matrix & A, const gene_matrix & B, gene_matrix & X, int N){
00088     X.noalias() = A*B;
00089   }
00090 
00091   static inline void transposed_matrix_matrix_product(const gene_matrix & A, const gene_matrix & B, gene_matrix & X, int N){
00092     X.noalias() = A.transpose()*B.transpose();
00093   }
00094 
00095   static inline void ata_product(const gene_matrix & A, gene_matrix & X, int N){
00096     X.noalias() = A.transpose()*A;
00097   }
00098 
00099   static inline void aat_product(const gene_matrix & A, gene_matrix & X, int N){
00100     X.noalias() = A*A.transpose();
00101   }
00102 
00103   static inline void matrix_vector_product(const gene_matrix & A, const gene_vector & B, gene_vector & X, int N){
00104     X.noalias() = A*B;
00105   }
00106 
00107   static inline void symv(const gene_matrix & A, const gene_vector & B, gene_vector & X, int N){
00108     X.noalias() = (A.template selfadjointView<Lower>() * B);
00109 //     internal::product_selfadjoint_vector<real,0,LowerTriangularBit,false,false>(N,A.data(),N, B.data(), 1, X.data(), 1);
00110   }
00111 
00112   template<typename Dest, typename Src> static void triassign(Dest& dst, const Src& src)
00113   {
00114     typedef typename Dest::Scalar Scalar;
00115     typedef typename internal::packet_traits<Scalar>::type Packet;
00116     const int PacketSize = sizeof(Packet)/sizeof(Scalar);
00117     int size = dst.cols();
00118     for(int j=0; j<size; j+=1)
00119     {
00120 //       const int alignedEnd = alignedStart + ((innerSize-alignedStart) & ~packetAlignedMask);
00121       Scalar* A0 = dst.data() + j*dst.stride();
00122       int starti = j;
00123       int alignedEnd = starti;
00124       int alignedStart = (starti) + internal::first_aligned(&A0[starti], size-starti);
00125       alignedEnd = alignedStart + ((size-alignedStart)/(2*PacketSize))*(PacketSize*2);
00126 
00127       // do the non-vectorizable part of the assignment
00128       for (int index = starti; index<alignedStart ; ++index)
00129       {
00130         if(Dest::Flags&RowMajorBit)
00131           dst.copyCoeff(j, index, src);
00132         else
00133           dst.copyCoeff(index, j, src);
00134       }
00135 
00136       // do the vectorizable part of the assignment
00137       for (int index = alignedStart; index<alignedEnd; index+=PacketSize)
00138       {
00139         if(Dest::Flags&RowMajorBit)
00140           dst.template copyPacket<Src, Aligned, Unaligned>(j, index, src);
00141         else
00142           dst.template copyPacket<Src, Aligned, Unaligned>(index, j, src);
00143       }
00144 
00145       // do the non-vectorizable part of the assignment
00146       for (int index = alignedEnd; index<size; ++index)
00147       {
00148         if(Dest::Flags&RowMajorBit)
00149           dst.copyCoeff(j, index, src);
00150         else
00151           dst.copyCoeff(index, j, src);
00152       }
00153       //dst.col(j).tail(N-j) = src.col(j).tail(N-j);
00154     }
00155   }
00156 
00157   static EIGEN_DONT_INLINE void syr2(gene_matrix & A,  gene_vector & X, gene_vector & Y, int N){
00158     // internal::product_selfadjoint_rank2_update<real,0,LowerTriangularBit>(N,A.data(),N, X.data(), 1, Y.data(), 1, -1);
00159     for(int j=0; j<N; ++j)
00160       A.col(j).tail(N-j) += X[j] * Y.tail(N-j) + Y[j] * X.tail(N-j);
00161   }
00162 
00163   static EIGEN_DONT_INLINE void ger(gene_matrix & A,  gene_vector & X, gene_vector & Y, int N){
00164     for(int j=0; j<N; ++j)
00165       A.col(j) += X * Y[j];
00166   }
00167 
00168   static EIGEN_DONT_INLINE void rot(gene_vector & A,  gene_vector & B, real c, real s, int N){
00169     internal::apply_rotation_in_the_plane(A, B, JacobiRotation<real>(c,s));
00170   }
00171 
00172   static inline void atv_product(gene_matrix & A, gene_vector & B, gene_vector & X, int N){
00173     X.noalias() = (A.transpose()*B);
00174   }
00175 
00176   static inline void axpy(real coef, const gene_vector & X, gene_vector & Y, int N){
00177     Y += coef * X;
00178   }
00179 
00180   static inline void axpby(real a, const gene_vector & X, real b, gene_vector & Y, int N){
00181     Y = a*X + b*Y;
00182   }
00183 
00184   static EIGEN_DONT_INLINE void copy_matrix(const gene_matrix & source, gene_matrix & cible, int N){
00185     cible = source;
00186   }
00187 
00188   static EIGEN_DONT_INLINE void copy_vector(const gene_vector & source, gene_vector & cible, int N){
00189     cible = source;
00190   }
00191 
00192   static inline void trisolve_lower(const gene_matrix & L, const gene_vector& B, gene_vector& X, int N){
00193     X = L.template triangularView<Lower>().solve(B);
00194   }
00195 
00196   static inline void trisolve_lower_matrix(const gene_matrix & L, const gene_matrix& B, gene_matrix& X, int N){
00197     X = L.template triangularView<Lower>().solve(B);
00198   }
00199 
00200   static inline void trmm(const gene_matrix & L, const gene_matrix& B, gene_matrix& X, int N){
00201     X = L.template triangularView<Lower>() * B;
00202   }
00203 
00204   static inline void cholesky(const gene_matrix & X, gene_matrix & C, int N){
00205     C = X;
00206     internal::llt_inplace<Lower>::blocked(C);
00207     //C = X.llt().matrixL();
00208 //     C = X;
00209 //     Cholesky<gene_matrix>::computeInPlace(C);
00210 //     Cholesky<gene_matrix>::computeInPlaceBlock(C);
00211   }
00212 
00213   static inline void lu_decomp(const gene_matrix & X, gene_matrix & C, int N){
00214     C = X.fullPivLu().matrixLU();
00215   }
00216 
00217   static inline void partial_lu_decomp(const gene_matrix & X, gene_matrix & C, int N){
00218     Matrix<DenseIndex,1,Dynamic> piv(N);
00219     DenseIndex nb;
00220     C = X;
00221     internal::partial_lu_inplace(C,piv,nb);
00222 //     C = X.partialPivLu().matrixLU();
00223   }
00224 
00225   static inline void tridiagonalization(const gene_matrix & X, gene_matrix & C, int N){
00226     typename Tridiagonalization<gene_matrix>::CoeffVectorType aux(N-1);
00227     C = X;
00228     internal::tridiagonalization_inplace(C, aux);
00229   }
00230 
00231   static inline void hessenberg(const gene_matrix & X, gene_matrix & C, int N){
00232     C = HessenbergDecomposition<gene_matrix>(X).packedMatrix();
00233   }
00234 
00235 
00236 
00237 };
00238 
00239 #endif


libicr
Author(s): Robert Krug
autogenerated on Mon Jan 6 2014 11:32:39