aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/nestbyvalue.cpp
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2019-01-17 18:27:25 +0100
committerGravatar Gael Guennebaud <g.gael@free.fr>2019-01-17 18:27:25 +0100
commit0fe6b7d687430fd1fe2d390da2f09fcb8ddc8093 (patch)
tree4d0aeef4a8a573f05bb61c284e5a48cca16bd473 /test/nestbyvalue.cpp
parent4b7cf7ff82a5bfa252dd2e00b449073272482d65 (diff)
Make nestByValue works again (broken since 3.3) and add unit tests.
Diffstat (limited to 'test/nestbyvalue.cpp')
-rw-r--r--test/nestbyvalue.cpp37
1 files changed, 37 insertions, 0 deletions
diff --git a/test/nestbyvalue.cpp b/test/nestbyvalue.cpp
new file mode 100644
index 000000000..c5356bc24
--- /dev/null
+++ b/test/nestbyvalue.cpp
@@ -0,0 +1,37 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2019 Gael Guennebaud <gael.guennebaud@inria.fr>
+//
+// This Source Code Form is subject to the terms of the Mozilla
+// Public License v. 2.0. If a copy of the MPL was not distributed
+// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
+
+#define TEST_ENABLE_TEMPORARY_TRACKING
+
+#include "main.h"
+
+typedef NestByValue<MatrixXd> CpyMatrixXd;
+typedef CwiseBinaryOp<internal::scalar_sum_op<double,double>,const CpyMatrixXd,const CpyMatrixXd> XprType;
+
+XprType get_xpr_with_temps(const MatrixXd& a)
+{
+ MatrixXd t1 = a.rowwise().reverse();
+ MatrixXd t2 = a+a;
+ return t1.nestByValue() + t2.nestByValue();
+}
+
+EIGEN_DECLARE_TEST(nestbyvalue)
+{
+ for(int i = 0; i < g_repeat; i++) {
+ Index rows = internal::random<Index>(1,EIGEN_TEST_MAX_SIZE);
+ Index cols = internal::random<Index>(1,EIGEN_TEST_MAX_SIZE);
+ MatrixXd a = MatrixXd(rows,cols);
+ nb_temporaries = 0;
+ XprType x = get_xpr_with_temps(a);
+ VERIFY_IS_EQUAL(nb_temporaries,6);
+ MatrixXd b = x;
+ VERIFY_IS_EQUAL(nb_temporaries,6+1);
+ VERIFY_IS_APPROX(b, a.rowwise().reverse().eval() + (a+a).eval());
+ }
+}