You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

240 lines
7.7 KiB

  1. //=====================================================
  2. // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
  3. //=====================================================
  4. //
  5. // This program is free software; you can redistribute it and/or
  6. // modify it under the terms of the GNU General Public License
  7. // as published by the Free Software Foundation; either version 2
  8. // of the License, or (at your option) any later version.
  9. //
  10. // This program is distributed in the hope that it will be useful,
  11. // but WITHOUT ANY WARRANTY; without even the implied warranty of
  12. // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  13. // GNU General Public License for more details.
  14. // You should have received a copy of the GNU General Public License
  15. // along with this program; if not, write to the Free Software
  16. // Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
  17. //
  18. #ifndef EIGEN3_INTERFACE_HH
  19. #define EIGEN3_INTERFACE_HH
  20. #include <Eigen/Eigen>
  21. #include <vector>
  22. #include "btl.hh"
  23. using namespace Eigen;
  24. template<class real, int SIZE=Dynamic>
  25. class eigen3_interface
  26. {
  27. public :
  28. enum {IsFixedSize = (SIZE!=Dynamic)};
  29. typedef real real_type;
  30. typedef std::vector<real> stl_vector;
  31. typedef std::vector<stl_vector> stl_matrix;
  32. typedef Eigen::Matrix<real,SIZE,SIZE> gene_matrix;
  33. typedef Eigen::Matrix<real,SIZE,1> gene_vector;
  34. static inline std::string name( void )
  35. {
  36. return EIGEN_MAKESTRING(BTL_PREFIX);
  37. }
  38. static void free_matrix(gene_matrix & A, int N) {}
  39. static void free_vector(gene_vector & B) {}
  40. static BTL_DONT_INLINE void matrix_from_stl(gene_matrix & A, stl_matrix & A_stl){
  41. A.resize(A_stl[0].size(), A_stl.size());
  42. for (int j=0; j<A_stl.size() ; j++){
  43. for (int i=0; i<A_stl[j].size() ; i++){
  44. A.coeffRef(i,j) = A_stl[j][i];
  45. }
  46. }
  47. }
  48. static BTL_DONT_INLINE void vector_from_stl(gene_vector & B, stl_vector & B_stl){
  49. B.resize(B_stl.size(),1);
  50. for (int i=0; i<B_stl.size() ; i++){
  51. B.coeffRef(i) = B_stl[i];
  52. }
  53. }
  54. static BTL_DONT_INLINE void vector_to_stl(gene_vector & B, stl_vector & B_stl){
  55. for (int i=0; i<B_stl.size() ; i++){
  56. B_stl[i] = B.coeff(i);
  57. }
  58. }
  59. static BTL_DONT_INLINE void matrix_to_stl(gene_matrix & A, stl_matrix & A_stl){
  60. int N=A_stl.size();
  61. for (int j=0;j<N;j++){
  62. A_stl[j].resize(N);
  63. for (int i=0;i<N;i++){
  64. A_stl[j][i] = A.coeff(i,j);
  65. }
  66. }
  67. }
  68. static inline void matrix_matrix_product(const gene_matrix & A, const gene_matrix & B, gene_matrix & X, int N){
  69. X.noalias() = A*B;
  70. }
  71. static inline void transposed_matrix_matrix_product(const gene_matrix & A, const gene_matrix & B, gene_matrix & X, int N){
  72. X.noalias() = A.transpose()*B.transpose();
  73. }
  74. // static inline void ata_product(const gene_matrix & A, gene_matrix & X, int N){
  75. // X.noalias() = A.transpose()*A;
  76. // }
  77. static inline void aat_product(const gene_matrix & A, gene_matrix & X, int N){
  78. X.template triangularView<Lower>().setZero();
  79. X.template selfadjointView<Lower>().rankUpdate(A);
  80. }
  81. static inline void matrix_vector_product(const gene_matrix & A, const gene_vector & B, gene_vector & X, int N){
  82. X.noalias() = A*B;
  83. }
  84. static inline void symv(const gene_matrix & A, const gene_vector & B, gene_vector & X, int N){
  85. X.noalias() = (A.template selfadjointView<Lower>() * B);
  86. // internal::product_selfadjoint_vector<real,0,LowerTriangularBit,false,false>(N,A.data(),N, B.data(), 1, X.data(), 1);
  87. }
  88. template<typename Dest, typename Src> static void triassign(Dest& dst, const Src& src)
  89. {
  90. typedef typename Dest::Scalar Scalar;
  91. typedef typename internal::packet_traits<Scalar>::type Packet;
  92. const int PacketSize = sizeof(Packet)/sizeof(Scalar);
  93. int size = dst.cols();
  94. for(int j=0; j<size; j+=1)
  95. {
  96. // const int alignedEnd = alignedStart + ((innerSize-alignedStart) & ~packetAlignedMask);
  97. Scalar* A0 = dst.data() + j*dst.stride();
  98. int starti = j;
  99. int alignedEnd = starti;
  100. int alignedStart = (starti) + internal::first_aligned(&A0[starti], size-starti);
  101. alignedEnd = alignedStart + ((size-alignedStart)/(2*PacketSize))*(PacketSize*2);
  102. // do the non-vectorizable part of the assignment
  103. for (int index = starti; index<alignedStart ; ++index)
  104. {
  105. if(Dest::Flags&RowMajorBit)
  106. dst.copyCoeff(j, index, src);
  107. else
  108. dst.copyCoeff(index, j, src);
  109. }
  110. // do the vectorizable part of the assignment
  111. for (int index = alignedStart; index<alignedEnd; index+=PacketSize)
  112. {
  113. if(Dest::Flags&RowMajorBit)
  114. dst.template copyPacket<Src, Aligned, Unaligned>(j, index, src);
  115. else
  116. dst.template copyPacket<Src, Aligned, Unaligned>(index, j, src);
  117. }
  118. // do the non-vectorizable part of the assignment
  119. for (int index = alignedEnd; index<size; ++index)
  120. {
  121. if(Dest::Flags&RowMajorBit)
  122. dst.copyCoeff(j, index, src);
  123. else
  124. dst.copyCoeff(index, j, src);
  125. }
  126. //dst.col(j).tail(N-j) = src.col(j).tail(N-j);
  127. }
  128. }
  129. static EIGEN_DONT_INLINE void syr2(gene_matrix & A, gene_vector & X, gene_vector & Y, int N){
  130. // internal::product_selfadjoint_rank2_update<real,0,LowerTriangularBit>(N,A.data(),N, X.data(), 1, Y.data(), 1, -1);
  131. for(int j=0; j<N; ++j)
  132. A.col(j).tail(N-j) += X[j] * Y.tail(N-j) + Y[j] * X.tail(N-j);
  133. }
  134. static EIGEN_DONT_INLINE void ger(gene_matrix & A, gene_vector & X, gene_vector & Y, int N){
  135. for(int j=0; j<N; ++j)
  136. A.col(j) += X * Y[j];
  137. }
  138. static EIGEN_DONT_INLINE void rot(gene_vector & A, gene_vector & B, real c, real s, int N){
  139. internal::apply_rotation_in_the_plane(A, B, JacobiRotation<real>(c,s));
  140. }
  141. static inline void atv_product(gene_matrix & A, gene_vector & B, gene_vector & X, int N){
  142. X.noalias() = (A.transpose()*B);
  143. }
  144. static inline void axpy(real coef, const gene_vector & X, gene_vector & Y, int N){
  145. Y += coef * X;
  146. }
  147. static inline void axpby(real a, const gene_vector & X, real b, gene_vector & Y, int N){
  148. Y = a*X + b*Y;
  149. }
  150. static EIGEN_DONT_INLINE void copy_matrix(const gene_matrix & source, gene_matrix & cible, int N){
  151. cible = source;
  152. }
  153. static EIGEN_DONT_INLINE void copy_vector(const gene_vector & source, gene_vector & cible, int N){
  154. cible = source;
  155. }
  156. static inline void trisolve_lower(const gene_matrix & L, const gene_vector& B, gene_vector& X, int N){
  157. X = L.template triangularView<Lower>().solve(B);
  158. }
  159. static inline void trisolve_lower_matrix(const gene_matrix & L, const gene_matrix& B, gene_matrix& X, int N){
  160. X = L.template triangularView<Upper>().solve(B);
  161. }
  162. static inline void trmm(const gene_matrix & L, const gene_matrix& B, gene_matrix& X, int N){
  163. X.noalias() = L.template triangularView<Lower>() * B;
  164. }
  165. static inline void cholesky(const gene_matrix & X, gene_matrix & C, int N){
  166. C = X;
  167. internal::llt_inplace<real,Lower>::blocked(C);
  168. //C = X.llt().matrixL();
  169. // C = X;
  170. // Cholesky<gene_matrix>::computeInPlace(C);
  171. // Cholesky<gene_matrix>::computeInPlaceBlock(C);
  172. }
  173. static inline void lu_decomp(const gene_matrix & X, gene_matrix & C, int N){
  174. C = X.fullPivLu().matrixLU();
  175. }
  176. static inline void partial_lu_decomp(const gene_matrix & X, gene_matrix & C, int N){
  177. Matrix<DenseIndex,1,Dynamic> piv(N);
  178. DenseIndex nb;
  179. C = X;
  180. internal::partial_lu_inplace(C,piv,nb);
  181. // C = X.partialPivLu().matrixLU();
  182. }
  183. static inline void tridiagonalization(const gene_matrix & X, gene_matrix & C, int N){
  184. typename Tridiagonalization<gene_matrix>::CoeffVectorType aux(N-1);
  185. C = X;
  186. internal::tridiagonalization_inplace(C, aux);
  187. }
  188. static inline void hessenberg(const gene_matrix & X, gene_matrix & C, int N){
  189. C = HessenbergDecomposition<gene_matrix>(X).packedMatrix();
  190. }
  191. };
  192. #endif