1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#define EIGEN_DEBUG_ASSIGN
11#include "main.h"
12#include <typeinfo>
13
14std::string demangle_traversal(int t)
15{
16  if(t==DefaultTraversal) return "DefaultTraversal";
17  if(t==LinearTraversal) return "LinearTraversal";
18  if(t==InnerVectorizedTraversal) return "InnerVectorizedTraversal";
19  if(t==LinearVectorizedTraversal) return "LinearVectorizedTraversal";
20  if(t==SliceVectorizedTraversal) return "SliceVectorizedTraversal";
21  return "?";
22}
23std::string demangle_unrolling(int t)
24{
25  if(t==NoUnrolling) return "NoUnrolling";
26  if(t==InnerUnrolling) return "InnerUnrolling";
27  if(t==CompleteUnrolling) return "CompleteUnrolling";
28  return "?";
29}
30
31template<typename Dst, typename Src>
32bool test_assign(const Dst&, const Src&, int traversal, int unrolling)
33{
34  internal::assign_traits<Dst,Src>::debug();
35  bool res = internal::assign_traits<Dst,Src>::Traversal==traversal
36          && internal::assign_traits<Dst,Src>::Unrolling==unrolling;
37  if(!res)
38  {
39    std::cerr << " Expected Traversal == " << demangle_traversal(traversal)
40              << " got " << demangle_traversal(internal::assign_traits<Dst,Src>::Traversal) << "\n";
41    std::cerr << " Expected Unrolling == " << demangle_unrolling(unrolling)
42              << " got " << demangle_unrolling(internal::assign_traits<Dst,Src>::Unrolling) << "\n";
43  }
44  return res;
45}
46
47template<typename Dst, typename Src>
48bool test_assign(int traversal, int unrolling)
49{
50  internal::assign_traits<Dst,Src>::debug();
51  bool res = internal::assign_traits<Dst,Src>::Traversal==traversal
52          && internal::assign_traits<Dst,Src>::Unrolling==unrolling;
53  if(!res)
54  {
55    std::cerr << " Expected Traversal == " << demangle_traversal(traversal)
56              << " got " << demangle_traversal(internal::assign_traits<Dst,Src>::Traversal) << "\n";
57    std::cerr << " Expected Unrolling == " << demangle_unrolling(unrolling)
58              << " got " << demangle_unrolling(internal::assign_traits<Dst,Src>::Unrolling) << "\n";
59  }
60  return res;
61}
62
63template<typename Xpr>
64bool test_redux(const Xpr&, int traversal, int unrolling)
65{
66  typedef internal::redux_traits<internal::scalar_sum_op<typename Xpr::Scalar>,Xpr> traits;
67  bool res = traits::Traversal==traversal && traits::Unrolling==unrolling;
68  if(!res)
69  {
70    std::cerr << " Expected Traversal == " << demangle_traversal(traversal)
71              << " got " << demangle_traversal(traits::Traversal) << "\n";
72    std::cerr << " Expected Unrolling == " << demangle_unrolling(unrolling)
73              << " got " << demangle_unrolling(traits::Unrolling) << "\n";
74  }
75  return res;
76}
77
78template<typename Scalar, bool Enable = internal::packet_traits<Scalar>::Vectorizable> struct vectorization_logic
79{
80  enum {
81    PacketSize = internal::packet_traits<Scalar>::size
82  };
83  static void run()
84  {
85
86    typedef Matrix<Scalar,PacketSize,1> Vector1;
87    typedef Matrix<Scalar,Dynamic,1> VectorX;
88    typedef Matrix<Scalar,Dynamic,Dynamic> MatrixXX;
89    typedef Matrix<Scalar,PacketSize,PacketSize> Matrix11;
90    typedef Matrix<Scalar,2*PacketSize,2*PacketSize> Matrix22;
91    typedef Matrix<Scalar,(Matrix11::Flags&RowMajorBit)?16:4*PacketSize,(Matrix11::Flags&RowMajorBit)?4*PacketSize:16> Matrix44;
92    typedef Matrix<Scalar,(Matrix11::Flags&RowMajorBit)?16:4*PacketSize,(Matrix11::Flags&RowMajorBit)?4*PacketSize:16,DontAlign|EIGEN_DEFAULT_MATRIX_STORAGE_ORDER_OPTION> Matrix44u;
93    typedef Matrix<Scalar,4*PacketSize,16,ColMajor> Matrix44c;
94    typedef Matrix<Scalar,4*PacketSize,16,RowMajor> Matrix44r;
95
96    typedef Matrix<Scalar,
97        (PacketSize==8 ? 4 : PacketSize==4 ? 2 : PacketSize==2 ? 1 : /*PacketSize==1 ?*/ 1),
98        (PacketSize==8 ? 2 : PacketSize==4 ? 2 : PacketSize==2 ? 2 : /*PacketSize==1 ?*/ 1)
99      > Matrix1;
100
101    typedef Matrix<Scalar,
102        (PacketSize==8 ? 4 : PacketSize==4 ? 2 : PacketSize==2 ? 1 : /*PacketSize==1 ?*/ 1),
103        (PacketSize==8 ? 2 : PacketSize==4 ? 2 : PacketSize==2 ? 2 : /*PacketSize==1 ?*/ 1),
104      DontAlign|((Matrix1::Flags&RowMajorBit)?RowMajor:ColMajor)> Matrix1u;
105
106    // this type is made such that it can only be vectorized when viewed as a linear 1D vector
107    typedef Matrix<Scalar,
108        (PacketSize==8 ? 4 : PacketSize==4 ? 6 : PacketSize==2 ? ((Matrix11::Flags&RowMajorBit)?2:3) : /*PacketSize==1 ?*/ 1),
109        (PacketSize==8 ? 6 : PacketSize==4 ? 2 : PacketSize==2 ? ((Matrix11::Flags&RowMajorBit)?3:2) : /*PacketSize==1 ?*/ 3)
110      > Matrix3;
111
112    #if !EIGEN_GCC_AND_ARCH_DOESNT_WANT_STACK_ALIGNMENT
113    VERIFY(test_assign(Vector1(),Vector1(),
114      InnerVectorizedTraversal,CompleteUnrolling));
115    VERIFY(test_assign(Vector1(),Vector1()+Vector1(),
116      InnerVectorizedTraversal,CompleteUnrolling));
117    VERIFY(test_assign(Vector1(),Vector1().cwiseProduct(Vector1()),
118      InnerVectorizedTraversal,CompleteUnrolling));
119    VERIFY(test_assign(Vector1(),Vector1().template cast<Scalar>(),
120      InnerVectorizedTraversal,CompleteUnrolling));
121
122
123    VERIFY(test_assign(Vector1(),Vector1(),
124      InnerVectorizedTraversal,CompleteUnrolling));
125    VERIFY(test_assign(Vector1(),Vector1()+Vector1(),
126      InnerVectorizedTraversal,CompleteUnrolling));
127    VERIFY(test_assign(Vector1(),Vector1().cwiseProduct(Vector1()),
128      InnerVectorizedTraversal,CompleteUnrolling));
129
130    VERIFY(test_assign(Matrix44(),Matrix44()+Matrix44(),
131      InnerVectorizedTraversal,InnerUnrolling));
132
133    VERIFY(test_assign(Matrix44u(),Matrix44()+Matrix44(),
134      LinearTraversal,NoUnrolling));
135
136    VERIFY(test_assign(Matrix1u(),Matrix1()+Matrix1(),
137      LinearTraversal,CompleteUnrolling));
138
139    VERIFY(test_assign(Matrix44c().col(1),Matrix44c().col(2)+Matrix44c().col(3),
140      InnerVectorizedTraversal,CompleteUnrolling));
141
142    VERIFY(test_assign(Matrix44r().row(2),Matrix44r().row(1)+Matrix44r().row(1),
143      InnerVectorizedTraversal,CompleteUnrolling));
144
145    if(PacketSize>1)
146    {
147      typedef Matrix<Scalar,3,3,ColMajor> Matrix33c;
148      VERIFY(test_assign(Matrix33c().row(2),Matrix33c().row(1)+Matrix33c().row(1),
149        LinearTraversal,CompleteUnrolling));
150      VERIFY(test_assign(Matrix33c().col(0),Matrix33c().col(1)+Matrix33c().col(1),
151        LinearTraversal,CompleteUnrolling));
152
153      VERIFY(test_assign(Matrix3(),Matrix3().cwiseQuotient(Matrix3()),
154        LinearVectorizedTraversal,CompleteUnrolling));
155
156      VERIFY(test_assign(Matrix<Scalar,17,17>(),Matrix<Scalar,17,17>()+Matrix<Scalar,17,17>(),
157        LinearTraversal,NoUnrolling));
158
159      VERIFY(test_assign(Matrix11(),Matrix<Scalar,17,17>().template block<PacketSize,PacketSize>(2,3)+Matrix<Scalar,17,17>().template block<PacketSize,PacketSize>(10,4),
160      DefaultTraversal,CompleteUnrolling));
161    }
162
163    VERIFY(test_redux(Matrix3(),
164      LinearVectorizedTraversal,CompleteUnrolling));
165
166    VERIFY(test_redux(Matrix44(),
167      LinearVectorizedTraversal,NoUnrolling));
168
169    VERIFY(test_redux(Matrix44().template block<(Matrix1::Flags&RowMajorBit)?4:PacketSize,(Matrix1::Flags&RowMajorBit)?PacketSize:4>(1,2),
170      DefaultTraversal,CompleteUnrolling));
171
172    VERIFY(test_redux(Matrix44c().template block<2*PacketSize,1>(1,2),
173      LinearVectorizedTraversal,CompleteUnrolling));
174
175    VERIFY(test_redux(Matrix44r().template block<1,2*PacketSize>(2,1),
176      LinearVectorizedTraversal,CompleteUnrolling));
177
178    VERIFY((test_assign<
179            Map<Matrix22, Aligned, OuterStride<3*PacketSize> >,
180            Matrix22
181            >(InnerVectorizedTraversal,CompleteUnrolling)));
182
183    VERIFY((test_assign<
184            Map<Matrix22, Aligned, InnerStride<3*PacketSize> >,
185            Matrix22
186            >(DefaultTraversal,CompleteUnrolling)));
187
188    VERIFY((test_assign(Matrix11(), Matrix11()*Matrix11(), InnerVectorizedTraversal, CompleteUnrolling)));
189    #endif
190
191    VERIFY(test_assign(MatrixXX(10,10),MatrixXX(20,20).block(10,10,2,3),
192      SliceVectorizedTraversal,NoUnrolling));
193
194    VERIFY(test_redux(VectorX(10),
195      LinearVectorizedTraversal,NoUnrolling));
196
197
198  }
199};
200
201template<typename Scalar> struct vectorization_logic<Scalar,false>
202{
203  static void run() {}
204};
205
206void test_vectorization_logic()
207{
208
209#ifdef EIGEN_VECTORIZE
210
211  CALL_SUBTEST( vectorization_logic<float>::run() );
212  CALL_SUBTEST( vectorization_logic<double>::run() );
213  CALL_SUBTEST( vectorization_logic<std::complex<float> >::run() );
214  CALL_SUBTEST( vectorization_logic<std::complex<double> >::run() );
215
216  if(internal::packet_traits<float>::Vectorizable)
217  {
218    VERIFY(test_assign(Matrix<float,3,3>(),Matrix<float,3,3>()+Matrix<float,3,3>(),
219      LinearTraversal,CompleteUnrolling));
220
221    VERIFY(test_redux(Matrix<float,5,2>(),
222      DefaultTraversal,CompleteUnrolling));
223  }
224
225  if(internal::packet_traits<double>::Vectorizable)
226  {
227    VERIFY(test_assign(Matrix<double,3,3>(),Matrix<double,3,3>()+Matrix<double,3,3>(),
228      LinearTraversal,CompleteUnrolling));
229
230    VERIFY(test_redux(Matrix<double,7,3>(),
231      DefaultTraversal,CompleteUnrolling));
232  }
233#endif // EIGEN_VECTORIZE
234
235}
236