Home | History | Annotate | Download | only in test
      1 
      2 #include "main.h"
      3 
      4 namespace Eigen {
      5 
      6   template<typename Lhs,typename Rhs>
      7   const Product<Lhs,Rhs>
      8   prod(const Lhs& lhs, const Rhs& rhs)
      9   {
     10     return Product<Lhs,Rhs>(lhs,rhs);
     11   }
     12 
     13   template<typename Lhs,typename Rhs>
     14   const Product<Lhs,Rhs,LazyProduct>
     15   lazyprod(const Lhs& lhs, const Rhs& rhs)
     16   {
     17     return Product<Lhs,Rhs,LazyProduct>(lhs,rhs);
     18   }
     19 
     20   template<typename DstXprType, typename SrcXprType>
     21   EIGEN_STRONG_INLINE
     22   DstXprType& copy_using_evaluator(const EigenBase<DstXprType> &dst, const SrcXprType &src)
     23   {
     24     call_assignment(dst.const_cast_derived(), src.derived(), internal::assign_op<typename DstXprType::Scalar,typename SrcXprType::Scalar>());
     25     return dst.const_cast_derived();
     26   }
     27 
     28   template<typename DstXprType, template <typename> class StorageBase, typename SrcXprType>
     29   EIGEN_STRONG_INLINE
     30   const DstXprType& copy_using_evaluator(const NoAlias<DstXprType, StorageBase>& dst, const SrcXprType &src)
     31   {
     32     call_assignment(dst, src.derived(), internal::assign_op<typename DstXprType::Scalar,typename SrcXprType::Scalar>());
     33     return dst.expression();
     34   }
     35 
     36   template<typename DstXprType, typename SrcXprType>
     37   EIGEN_STRONG_INLINE
     38   DstXprType& copy_using_evaluator(const PlainObjectBase<DstXprType> &dst, const SrcXprType &src)
     39   {
     40     #ifdef EIGEN_NO_AUTOMATIC_RESIZING
     41     eigen_assert((dst.size()==0 || (IsVectorAtCompileTime ? (dst.size() == src.size())
     42                                                           : (dst.rows() == src.rows() && dst.cols() == src.cols())))
     43                 && "Size mismatch. Automatic resizing is disabled because EIGEN_NO_AUTOMATIC_RESIZING is defined");
     44   #else
     45     dst.const_cast_derived().resizeLike(src.derived());
     46   #endif
     47 
     48     call_assignment(dst.const_cast_derived(), src.derived(), internal::assign_op<typename DstXprType::Scalar,typename SrcXprType::Scalar>());
     49     return dst.const_cast_derived();
     50   }
     51 
     52   template<typename DstXprType, typename SrcXprType>
     53   void add_assign_using_evaluator(const DstXprType& dst, const SrcXprType& src)
     54   {
     55     typedef typename DstXprType::Scalar Scalar;
     56     call_assignment(const_cast<DstXprType&>(dst), src.derived(), internal::add_assign_op<Scalar,typename SrcXprType::Scalar>());
     57   }
     58 
     59   template<typename DstXprType, typename SrcXprType>
     60   void subtract_assign_using_evaluator(const DstXprType& dst, const SrcXprType& src)
     61   {
     62     typedef typename DstXprType::Scalar Scalar;
     63     call_assignment(const_cast<DstXprType&>(dst), src.derived(), internal::sub_assign_op<Scalar,typename SrcXprType::Scalar>());
     64   }
     65 
     66   template<typename DstXprType, typename SrcXprType>
     67   void multiply_assign_using_evaluator(const DstXprType& dst, const SrcXprType& src)
     68   {
     69     typedef typename DstXprType::Scalar Scalar;
     70     call_assignment(dst.const_cast_derived(), src.derived(), internal::mul_assign_op<Scalar,typename SrcXprType::Scalar>());
     71   }
     72 
     73   template<typename DstXprType, typename SrcXprType>
     74   void divide_assign_using_evaluator(const DstXprType& dst, const SrcXprType& src)
     75   {
     76     typedef typename DstXprType::Scalar Scalar;
     77     call_assignment(dst.const_cast_derived(), src.derived(), internal::div_assign_op<Scalar,typename SrcXprType::Scalar>());
     78   }
     79 
     80   template<typename DstXprType, typename SrcXprType>
     81   void swap_using_evaluator(const DstXprType& dst, const SrcXprType& src)
     82   {
     83     typedef typename DstXprType::Scalar Scalar;
     84     call_assignment(dst.const_cast_derived(), src.const_cast_derived(), internal::swap_assign_op<Scalar>());
     85   }
     86 
     87   namespace internal {
     88     template<typename Dst, template <typename> class StorageBase, typename Src, typename Func>
     89     EIGEN_DEVICE_FUNC void call_assignment(const NoAlias<Dst,StorageBase>& dst, const Src& src, const Func& func)
     90     {
     91       call_assignment_no_alias(dst.expression(), src, func);
     92     }
     93   }
     94 
     95 }
     96 
     97 template<typename XprType> long get_cost(const XprType& ) { return Eigen::internal::evaluator<XprType>::CoeffReadCost; }
     98 
     99 using namespace std;
    100 
    101 #define VERIFY_IS_APPROX_EVALUATOR(DEST,EXPR) VERIFY_IS_APPROX(copy_using_evaluator(DEST,(EXPR)), (EXPR).eval());
    102 #define VERIFY_IS_APPROX_EVALUATOR2(DEST,EXPR,REF) VERIFY_IS_APPROX(copy_using_evaluator(DEST,(EXPR)), (REF).eval());
    103 
    104 void test_evaluators()
    105 {
    106   // Testing Matrix evaluator and Transpose
    107   Vector2d v = Vector2d::Random();
    108   const Vector2d v_const(v);
    109   Vector2d v2;
    110   RowVector2d w;
    111 
    112   VERIFY_IS_APPROX_EVALUATOR(v2, v);
    113   VERIFY_IS_APPROX_EVALUATOR(v2, v_const);
    114 
    115   // Testing Transpose
    116   VERIFY_IS_APPROX_EVALUATOR(w, v.transpose()); // Transpose as rvalue
    117   VERIFY_IS_APPROX_EVALUATOR(w, v_const.transpose());
    118 
    119   copy_using_evaluator(w.transpose(), v); // Transpose as lvalue
    120   VERIFY_IS_APPROX(w,v.transpose().eval());
    121 
    122   copy_using_evaluator(w.transpose(), v_const);
    123   VERIFY_IS_APPROX(w,v_const.transpose().eval());
    124 
    125   // Testing Array evaluator
    126   {
    127     ArrayXXf a(2,3);
    128     ArrayXXf b(3,2);
    129     a << 1,2,3, 4,5,6;
    130     const ArrayXXf a_const(a);
    131 
    132     VERIFY_IS_APPROX_EVALUATOR(b, a.transpose());
    133 
    134     VERIFY_IS_APPROX_EVALUATOR(b, a_const.transpose());
    135 
    136     // Testing CwiseNullaryOp evaluator
    137     copy_using_evaluator(w, RowVector2d::Random());
    138     VERIFY((w.array() >= -1).all() && (w.array() <= 1).all()); // not easy to test ...
    139 
    140     VERIFY_IS_APPROX_EVALUATOR(w, RowVector2d::Zero());
    141 
    142     VERIFY_IS_APPROX_EVALUATOR(w, RowVector2d::Constant(3));
    143 
    144     // mix CwiseNullaryOp and transpose
    145     VERIFY_IS_APPROX_EVALUATOR(w, Vector2d::Zero().transpose());
    146   }
    147 
    148   {
    149     // test product expressions
    150     int s = internal::random<int>(1,100);
    151     MatrixXf a(s,s), b(s,s), c(s,s), d(s,s);
    152     a.setRandom();
    153     b.setRandom();
    154     c.setRandom();
    155     d.setRandom();
    156     VERIFY_IS_APPROX_EVALUATOR(d, (a + b));
    157     VERIFY_IS_APPROX_EVALUATOR(d, (a + b).transpose());
    158     VERIFY_IS_APPROX_EVALUATOR2(d, prod(a,b), a*b);
    159     VERIFY_IS_APPROX_EVALUATOR2(d.noalias(), prod(a,b), a*b);
    160     VERIFY_IS_APPROX_EVALUATOR2(d, prod(a,b) + c, a*b + c);
    161     VERIFY_IS_APPROX_EVALUATOR2(d, s * prod(a,b), s * a*b);
    162     VERIFY_IS_APPROX_EVALUATOR2(d, prod(a,b).transpose(), (a*b).transpose());
    163     VERIFY_IS_APPROX_EVALUATOR2(d, prod(a,b) + prod(b,c), a*b + b*c);
    164 
    165     // check that prod works even with aliasing present
    166     c = a*a;
    167     copy_using_evaluator(a, prod(a,a));
    168     VERIFY_IS_APPROX(a,c);
    169 
    170     // check compound assignment of products
    171     d = c;
    172     add_assign_using_evaluator(c.noalias(), prod(a,b));
    173     d.noalias() += a*b;
    174     VERIFY_IS_APPROX(c, d);
    175 
    176     d = c;
    177     subtract_assign_using_evaluator(c.noalias(), prod(a,b));
    178     d.noalias() -= a*b;
    179     VERIFY_IS_APPROX(c, d);
    180   }
    181 
    182   {
    183     // test product with all possible sizes
    184     int s = internal::random<int>(1,100);
    185     Matrix<float,      1,      1> m11, res11;  m11.setRandom(1,1);
    186     Matrix<float,      1,      4> m14, res14;  m14.setRandom(1,4);
    187     Matrix<float,      1,Dynamic> m1X, res1X;  m1X.setRandom(1,s);
    188     Matrix<float,      4,      1> m41, res41;  m41.setRandom(4,1);
    189     Matrix<float,      4,      4> m44, res44;  m44.setRandom(4,4);
    190     Matrix<float,      4,Dynamic> m4X, res4X;  m4X.setRandom(4,s);
    191     Matrix<float,Dynamic,      1> mX1, resX1;  mX1.setRandom(s,1);
    192     Matrix<float,Dynamic,      4> mX4, resX4;  mX4.setRandom(s,4);
    193     Matrix<float,Dynamic,Dynamic> mXX, resXX;  mXX.setRandom(s,s);
    194 
    195     VERIFY_IS_APPROX_EVALUATOR2(res11, prod(m11,m11), m11*m11);
    196     VERIFY_IS_APPROX_EVALUATOR2(res11, prod(m14,m41), m14*m41);
    197     VERIFY_IS_APPROX_EVALUATOR2(res11, prod(m1X,mX1), m1X*mX1);
    198     VERIFY_IS_APPROX_EVALUATOR2(res14, prod(m11,m14), m11*m14);
    199     VERIFY_IS_APPROX_EVALUATOR2(res14, prod(m14,m44), m14*m44);
    200     VERIFY_IS_APPROX_EVALUATOR2(res14, prod(m1X,mX4), m1X*mX4);
    201     VERIFY_IS_APPROX_EVALUATOR2(res1X, prod(m11,m1X), m11*m1X);
    202     VERIFY_IS_APPROX_EVALUATOR2(res1X, prod(m14,m4X), m14*m4X);
    203     VERIFY_IS_APPROX_EVALUATOR2(res1X, prod(m1X,mXX), m1X*mXX);
    204     VERIFY_IS_APPROX_EVALUATOR2(res41, prod(m41,m11), m41*m11);
    205     VERIFY_IS_APPROX_EVALUATOR2(res41, prod(m44,m41), m44*m41);
    206     VERIFY_IS_APPROX_EVALUATOR2(res41, prod(m4X,mX1), m4X*mX1);
    207     VERIFY_IS_APPROX_EVALUATOR2(res44, prod(m41,m14), m41*m14);
    208     VERIFY_IS_APPROX_EVALUATOR2(res44, prod(m44,m44), m44*m44);
    209     VERIFY_IS_APPROX_EVALUATOR2(res44, prod(m4X,mX4), m4X*mX4);
    210     VERIFY_IS_APPROX_EVALUATOR2(res4X, prod(m41,m1X), m41*m1X);
    211     VERIFY_IS_APPROX_EVALUATOR2(res4X, prod(m44,m4X), m44*m4X);
    212     VERIFY_IS_APPROX_EVALUATOR2(res4X, prod(m4X,mXX), m4X*mXX);
    213     VERIFY_IS_APPROX_EVALUATOR2(resX1, prod(mX1,m11), mX1*m11);
    214     VERIFY_IS_APPROX_EVALUATOR2(resX1, prod(mX4,m41), mX4*m41);
    215     VERIFY_IS_APPROX_EVALUATOR2(resX1, prod(mXX,mX1), mXX*mX1);
    216     VERIFY_IS_APPROX_EVALUATOR2(resX4, prod(mX1,m14), mX1*m14);
    217     VERIFY_IS_APPROX_EVALUATOR2(resX4, prod(mX4,m44), mX4*m44);
    218     VERIFY_IS_APPROX_EVALUATOR2(resX4, prod(mXX,mX4), mXX*mX4);
    219     VERIFY_IS_APPROX_EVALUATOR2(resXX, prod(mX1,m1X), mX1*m1X);
    220     VERIFY_IS_APPROX_EVALUATOR2(resXX, prod(mX4,m4X), mX4*m4X);
    221     VERIFY_IS_APPROX_EVALUATOR2(resXX, prod(mXX,mXX), mXX*mXX);
    222   }
    223 
    224   {
    225     ArrayXXf a(2,3);
    226     ArrayXXf b(3,2);
    227     a << 1,2,3, 4,5,6;
    228     const ArrayXXf a_const(a);
    229 
    230     // this does not work because Random is eval-before-nested:
    231     // copy_using_evaluator(w, Vector2d::Random().transpose());
    232 
    233     // test CwiseUnaryOp
    234     VERIFY_IS_APPROX_EVALUATOR(v2, 3 * v);
    235     VERIFY_IS_APPROX_EVALUATOR(w, (3 * v).transpose());
    236     VERIFY_IS_APPROX_EVALUATOR(b, (a + 3).transpose());
    237     VERIFY_IS_APPROX_EVALUATOR(b, (2 * a_const + 3).transpose());
    238 
    239     // test CwiseBinaryOp
    240     VERIFY_IS_APPROX_EVALUATOR(v2, v + Vector2d::Ones());
    241     VERIFY_IS_APPROX_EVALUATOR(w, (v + Vector2d::Ones()).transpose().cwiseProduct(RowVector2d::Constant(3)));
    242 
    243     // dynamic matrices and arrays
    244     MatrixXd mat1(6,6), mat2(6,6);
    245     VERIFY_IS_APPROX_EVALUATOR(mat1, MatrixXd::Identity(6,6));
    246     VERIFY_IS_APPROX_EVALUATOR(mat2, mat1);
    247     copy_using_evaluator(mat2.transpose(), mat1);
    248     VERIFY_IS_APPROX(mat2.transpose(), mat1);
    249 
    250     ArrayXXd arr1(6,6), arr2(6,6);
    251     VERIFY_IS_APPROX_EVALUATOR(arr1, ArrayXXd::Constant(6,6, 3.0));
    252     VERIFY_IS_APPROX_EVALUATOR(arr2, arr1);
    253 
    254     // test automatic resizing
    255     mat2.resize(3,3);
    256     VERIFY_IS_APPROX_EVALUATOR(mat2, mat1);
    257     arr2.resize(9,9);
    258     VERIFY_IS_APPROX_EVALUATOR(arr2, arr1);
    259 
    260     // test direct traversal
    261     Matrix3f m3;
    262     Array33f a3;
    263     VERIFY_IS_APPROX_EVALUATOR(m3, Matrix3f::Identity());  // matrix, nullary
    264     // TODO: find a way to test direct traversal with array
    265     VERIFY_IS_APPROX_EVALUATOR(m3.transpose(), Matrix3f::Identity().transpose());  // transpose
    266     VERIFY_IS_APPROX_EVALUATOR(m3, 2 * Matrix3f::Identity());  // unary
    267     VERIFY_IS_APPROX_EVALUATOR(m3, Matrix3f::Identity() + Matrix3f::Zero());  // binary
    268     VERIFY_IS_APPROX_EVALUATOR(m3.block(0,0,2,2), Matrix3f::Identity().block(1,1,2,2));  // block
    269 
    270     // test linear traversal
    271     VERIFY_IS_APPROX_EVALUATOR(m3, Matrix3f::Zero());  // matrix, nullary
    272     VERIFY_IS_APPROX_EVALUATOR(a3, Array33f::Zero());  // array
    273     VERIFY_IS_APPROX_EVALUATOR(m3.transpose(), Matrix3f::Zero().transpose());  // transpose
    274     VERIFY_IS_APPROX_EVALUATOR(m3, 2 * Matrix3f::Zero());  // unary
    275     VERIFY_IS_APPROX_EVALUATOR(m3, Matrix3f::Zero() + m3);  // binary
    276 
    277     // test inner vectorization
    278     Matrix4f m4, m4src = Matrix4f::Random();
    279     Array44f a4, a4src = Matrix4f::Random();
    280     VERIFY_IS_APPROX_EVALUATOR(m4, m4src);  // matrix
    281     VERIFY_IS_APPROX_EVALUATOR(a4, a4src);  // array
    282     VERIFY_IS_APPROX_EVALUATOR(m4.transpose(), m4src.transpose());  // transpose
    283     // TODO: find out why Matrix4f::Zero() does not allow inner vectorization
    284     VERIFY_IS_APPROX_EVALUATOR(m4, 2 * m4src);  // unary
    285     VERIFY_IS_APPROX_EVALUATOR(m4, m4src + m4src);  // binary
    286 
    287     // test linear vectorization
    288     MatrixXf mX(6,6), mXsrc = MatrixXf::Random(6,6);
    289     ArrayXXf aX(6,6), aXsrc = ArrayXXf::Random(6,6);
    290     VERIFY_IS_APPROX_EVALUATOR(mX, mXsrc);  // matrix
    291     VERIFY_IS_APPROX_EVALUATOR(aX, aXsrc);  // array
    292     VERIFY_IS_APPROX_EVALUATOR(mX.transpose(), mXsrc.transpose());  // transpose
    293     VERIFY_IS_APPROX_EVALUATOR(mX, MatrixXf::Zero(6,6));  // nullary
    294     VERIFY_IS_APPROX_EVALUATOR(mX, 2 * mXsrc);  // unary
    295     VERIFY_IS_APPROX_EVALUATOR(mX, mXsrc + mXsrc);  // binary
    296 
    297     // test blocks and slice vectorization
    298     VERIFY_IS_APPROX_EVALUATOR(m4, (mXsrc.block<4,4>(1,0)));
    299     VERIFY_IS_APPROX_EVALUATOR(aX, ArrayXXf::Constant(10, 10, 3.0).block(2, 3, 6, 6));
    300 
    301     Matrix4f m4ref = m4;
    302     copy_using_evaluator(m4.block(1, 1, 2, 3), m3.bottomRows(2));
    303     m4ref.block(1, 1, 2, 3) = m3.bottomRows(2);
    304     VERIFY_IS_APPROX(m4, m4ref);
    305 
    306     mX.setIdentity(20,20);
    307     MatrixXf mXref = MatrixXf::Identity(20,20);
    308     mXsrc = MatrixXf::Random(9,12);
    309     copy_using_evaluator(mX.block(4, 4, 9, 12), mXsrc);
    310     mXref.block(4, 4, 9, 12) = mXsrc;
    311     VERIFY_IS_APPROX(mX, mXref);
    312 
    313     // test Map
    314     const float raw[3] = {1,2,3};
    315     float buffer[3] = {0,0,0};
    316     Vector3f v3;
    317     Array3f a3f;
    318     VERIFY_IS_APPROX_EVALUATOR(v3, Map<const Vector3f>(raw));
    319     VERIFY_IS_APPROX_EVALUATOR(a3f, Map<const Array3f>(raw));
    320     Vector3f::Map(buffer) = 2*v3;
    321     VERIFY(buffer[0] == 2);
    322     VERIFY(buffer[1] == 4);
    323     VERIFY(buffer[2] == 6);
    324 
    325     // test CwiseUnaryView
    326     mat1.setRandom();
    327     mat2.setIdentity();
    328     MatrixXcd matXcd(6,6), matXcd_ref(6,6);
    329     copy_using_evaluator(matXcd.real(), mat1);
    330     copy_using_evaluator(matXcd.imag(), mat2);
    331     matXcd_ref.real() = mat1;
    332     matXcd_ref.imag() = mat2;
    333     VERIFY_IS_APPROX(matXcd, matXcd_ref);
    334 
    335     // test Select
    336     VERIFY_IS_APPROX_EVALUATOR(aX, (aXsrc > 0).select(aXsrc, -aXsrc));
    337 
    338     // test Replicate
    339     mXsrc = MatrixXf::Random(6, 6);
    340     VectorXf vX = VectorXf::Random(6);
    341     mX.resize(6, 6);
    342     VERIFY_IS_APPROX_EVALUATOR(mX, mXsrc.colwise() + vX);
    343     matXcd.resize(12, 12);
    344     VERIFY_IS_APPROX_EVALUATOR(matXcd, matXcd_ref.replicate(2,2));
    345     VERIFY_IS_APPROX_EVALUATOR(matXcd, (matXcd_ref.replicate<2,2>()));
    346 
    347     // test partial reductions
    348     VectorXd vec1(6);
    349     VERIFY_IS_APPROX_EVALUATOR(vec1, mat1.rowwise().sum());
    350     VERIFY_IS_APPROX_EVALUATOR(vec1, mat1.colwise().sum().transpose());
    351 
    352     // test MatrixWrapper and ArrayWrapper
    353     mat1.setRandom(6,6);
    354     arr1.setRandom(6,6);
    355     VERIFY_IS_APPROX_EVALUATOR(mat2, arr1.matrix());
    356     VERIFY_IS_APPROX_EVALUATOR(arr2, mat1.array());
    357     VERIFY_IS_APPROX_EVALUATOR(mat2, (arr1 + 2).matrix());
    358     VERIFY_IS_APPROX_EVALUATOR(arr2, mat1.array() + 2);
    359     mat2.array() = arr1 * arr1;
    360     VERIFY_IS_APPROX(mat2, (arr1 * arr1).matrix());
    361     arr2.matrix() = MatrixXd::Identity(6,6);
    362     VERIFY_IS_APPROX(arr2, MatrixXd::Identity(6,6).array());
    363 
    364     // test Reverse
    365     VERIFY_IS_APPROX_EVALUATOR(arr2, arr1.reverse());
    366     VERIFY_IS_APPROX_EVALUATOR(arr2, arr1.colwise().reverse());
    367     VERIFY_IS_APPROX_EVALUATOR(arr2, arr1.rowwise().reverse());
    368     arr2.reverse() = arr1;
    369     VERIFY_IS_APPROX(arr2, arr1.reverse());
    370     mat2.array() = mat1.array().reverse();
    371     VERIFY_IS_APPROX(mat2.array(), mat1.array().reverse());
    372 
    373     // test Diagonal
    374     VERIFY_IS_APPROX_EVALUATOR(vec1, mat1.diagonal());
    375     vec1.resize(5);
    376     VERIFY_IS_APPROX_EVALUATOR(vec1, mat1.diagonal(1));
    377     VERIFY_IS_APPROX_EVALUATOR(vec1, mat1.diagonal<-1>());
    378     vec1.setRandom();
    379 
    380     mat2 = mat1;
    381     copy_using_evaluator(mat1.diagonal(1), vec1);
    382     mat2.diagonal(1) = vec1;
    383     VERIFY_IS_APPROX(mat1, mat2);
    384 
    385     copy_using_evaluator(mat1.diagonal<-1>(), mat1.diagonal(1));
    386     mat2.diagonal<-1>() = mat2.diagonal(1);
    387     VERIFY_IS_APPROX(mat1, mat2);
    388   }
    389 
    390   {
    391     // test swapping
    392     MatrixXd mat1, mat2, mat1ref, mat2ref;
    393     mat1ref = mat1 = MatrixXd::Random(6, 6);
    394     mat2ref = mat2 = 2 * mat1 + MatrixXd::Identity(6, 6);
    395     swap_using_evaluator(mat1, mat2);
    396     mat1ref.swap(mat2ref);
    397     VERIFY_IS_APPROX(mat1, mat1ref);
    398     VERIFY_IS_APPROX(mat2, mat2ref);
    399 
    400     swap_using_evaluator(mat1.block(0, 0, 3, 3), mat2.block(3, 3, 3, 3));
    401     mat1ref.block(0, 0, 3, 3).swap(mat2ref.block(3, 3, 3, 3));
    402     VERIFY_IS_APPROX(mat1, mat1ref);
    403     VERIFY_IS_APPROX(mat2, mat2ref);
    404 
    405     swap_using_evaluator(mat1.row(2), mat2.col(3).transpose());
    406     mat1.row(2).swap(mat2.col(3).transpose());
    407     VERIFY_IS_APPROX(mat1, mat1ref);
    408     VERIFY_IS_APPROX(mat2, mat2ref);
    409   }
    410 
    411   {
    412     // test compound assignment
    413     const Matrix4d mat_const = Matrix4d::Random();
    414     Matrix4d mat, mat_ref;
    415     mat = mat_ref = Matrix4d::Identity();
    416     add_assign_using_evaluator(mat, mat_const);
    417     mat_ref += mat_const;
    418     VERIFY_IS_APPROX(mat, mat_ref);
    419 
    420     subtract_assign_using_evaluator(mat.row(1), 2*mat.row(2));
    421     mat_ref.row(1) -= 2*mat_ref.row(2);
    422     VERIFY_IS_APPROX(mat, mat_ref);
    423 
    424     const ArrayXXf arr_const = ArrayXXf::Random(5,3);
    425     ArrayXXf arr, arr_ref;
    426     arr = arr_ref = ArrayXXf::Constant(5, 3, 0.5);
    427     multiply_assign_using_evaluator(arr, arr_const);
    428     arr_ref *= arr_const;
    429     VERIFY_IS_APPROX(arr, arr_ref);
    430 
    431     divide_assign_using_evaluator(arr.row(1), arr.row(2) + 1);
    432     arr_ref.row(1) /= (arr_ref.row(2) + 1);
    433     VERIFY_IS_APPROX(arr, arr_ref);
    434   }
    435 
    436   {
    437     // test triangular shapes
    438     MatrixXd A = MatrixXd::Random(6,6), B(6,6), C(6,6), D(6,6);
    439     A.setRandom();B.setRandom();
    440     VERIFY_IS_APPROX_EVALUATOR2(B, A.triangularView<Upper>(), MatrixXd(A.triangularView<Upper>()));
    441 
    442     A.setRandom();B.setRandom();
    443     VERIFY_IS_APPROX_EVALUATOR2(B, A.triangularView<UnitLower>(), MatrixXd(A.triangularView<UnitLower>()));
    444 
    445     A.setRandom();B.setRandom();
    446     VERIFY_IS_APPROX_EVALUATOR2(B, A.triangularView<UnitUpper>(), MatrixXd(A.triangularView<UnitUpper>()));
    447 
    448     A.setRandom();B.setRandom();
    449     C = B; C.triangularView<Upper>() = A;
    450     copy_using_evaluator(B.triangularView<Upper>(), A);
    451     VERIFY(B.isApprox(C) && "copy_using_evaluator(B.triangularView<Upper>(), A)");
    452 
    453     A.setRandom();B.setRandom();
    454     C = B; C.triangularView<Lower>() = A.triangularView<Lower>();
    455     copy_using_evaluator(B.triangularView<Lower>(), A.triangularView<Lower>());
    456     VERIFY(B.isApprox(C) && "copy_using_evaluator(B.triangularView<Lower>(), A.triangularView<Lower>())");
    457 
    458 
    459     A.setRandom();B.setRandom();
    460     C = B; C.triangularView<Lower>() = A.triangularView<Upper>().transpose();
    461     copy_using_evaluator(B.triangularView<Lower>(), A.triangularView<Upper>().transpose());
    462     VERIFY(B.isApprox(C) && "copy_using_evaluator(B.triangularView<Lower>(), A.triangularView<Lower>().transpose())");
    463 
    464 
    465     A.setRandom();B.setRandom(); C = B; D = A;
    466     C.triangularView<Upper>().swap(D.triangularView<Upper>());
    467     swap_using_evaluator(B.triangularView<Upper>(), A.triangularView<Upper>());
    468     VERIFY(B.isApprox(C) && "swap_using_evaluator(B.triangularView<Upper>(), A.triangularView<Upper>())");
    469 
    470 
    471     VERIFY_IS_APPROX_EVALUATOR2(B, prod(A.triangularView<Upper>(),A), MatrixXd(A.triangularView<Upper>()*A));
    472 
    473     VERIFY_IS_APPROX_EVALUATOR2(B, prod(A.selfadjointView<Upper>(),A), MatrixXd(A.selfadjointView<Upper>()*A));
    474   }
    475 
    476   {
    477     // test diagonal shapes
    478     VectorXd d = VectorXd::Random(6);
    479     MatrixXd A = MatrixXd::Random(6,6), B(6,6);
    480     A.setRandom();B.setRandom();
    481 
    482     VERIFY_IS_APPROX_EVALUATOR2(B, lazyprod(d.asDiagonal(),A), MatrixXd(d.asDiagonal()*A));
    483     VERIFY_IS_APPROX_EVALUATOR2(B, lazyprod(A,d.asDiagonal()), MatrixXd(A*d.asDiagonal()));
    484   }
    485 
    486   {
    487     // test CoeffReadCost
    488     Matrix4d a, b;
    489     VERIFY_IS_EQUAL( get_cost(a), 1 );
    490     VERIFY_IS_EQUAL( get_cost(a+b), 3);
    491     VERIFY_IS_EQUAL( get_cost(2*a+b), 4);
    492     VERIFY_IS_EQUAL( get_cost(a*b), 1);
    493     VERIFY_IS_EQUAL( get_cost(a.lazyProduct(b)), 15);
    494     VERIFY_IS_EQUAL( get_cost(a*(a*b)), 1);
    495     VERIFY_IS_EQUAL( get_cost(a.lazyProduct(a*b)), 15);
    496     VERIFY_IS_EQUAL( get_cost(a*(a+b)), 1);
    497     VERIFY_IS_EQUAL( get_cost(a.lazyProduct(a+b)), 15);
    498   }
    499 }
    500