-
Notifications
You must be signed in to change notification settings - Fork 477
/
Copy pathnesting_ops.cpp
121 lines (105 loc) · 4.53 KB
/
nesting_ops.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2010 Hauke Heibel <hauke.heibel@gmail.com>
// Copyright (C) 2015 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"
template <int N, typename XprType>
void use_n_times(const XprType& xpr) {
typename internal::nested_eval<XprType, N>::type mat(xpr);
typename XprType::PlainObject res(mat.rows(), mat.cols());
nb_temporaries--; // remove res
res.setZero();
for (int i = 0; i < N; ++i) res += mat;
}
template <int N, typename ReferenceType, typename XprType>
bool verify_eval_type(const XprType&, const ReferenceType&) {
typedef typename internal::nested_eval<XprType, N>::type EvalType;
return internal::is_same<
typename internal::remove_all<EvalType>::type,
typename internal::remove_all<ReferenceType>::type>::value;
}
template <typename MatrixType>
void run_nesting_ops_1(const MatrixType& _m) {
typename internal::nested_eval<MatrixType, 2>::type m(_m);
// Make really sure that we are in debug mode!
VERIFY_RAISES_ASSERT(eigen_assert(false));
// The only intention of these tests is to ensure that this code does
// not trigger any asserts or segmentation faults... more to come.
VERIFY_IS_APPROX((m.transpose() * m).diagonal().sum(),
(m.transpose() * m).diagonal().sum());
VERIFY_IS_APPROX((m.transpose() * m).diagonal().array().abs().sum(),
(m.transpose() * m).diagonal().array().abs().sum());
VERIFY_IS_APPROX((m.transpose() * m).array().abs().sum(),
(m.transpose() * m).array().abs().sum());
}
template <typename MatrixType>
void run_nesting_ops_2(const MatrixType& _m) {
typedef typename MatrixType::Scalar Scalar;
Index rows = _m.rows();
Index cols = _m.cols();
MatrixType m1 = MatrixType::Random(rows, cols);
Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime,
ColMajor>
m2;
if ((MatrixType::SizeAtCompileTime == Dynamic)) {
VERIFY_EVALUATION_COUNT(use_n_times<1>(m1 + m1 * m1), 1);
VERIFY_EVALUATION_COUNT(use_n_times<10>(m1 + m1 * m1), 1);
VERIFY_EVALUATION_COUNT(
use_n_times<1>(m1.template triangularView<Lower>().solve(m1.col(0))),
1);
VERIFY_EVALUATION_COUNT(
use_n_times<10>(m1.template triangularView<Lower>().solve(m1.col(0))),
1);
VERIFY_EVALUATION_COUNT(
use_n_times<1>(Scalar(2) *
m1.template triangularView<Lower>().solve(m1.col(0))),
2); // FIXME could be one by applying the scaling in-place on the solve
// result
VERIFY_EVALUATION_COUNT(
use_n_times<1>(m1.col(0) +
m1.template triangularView<Lower>().solve(m1.col(0))),
2); // FIXME could be one by adding m1.col() inplace
VERIFY_EVALUATION_COUNT(
use_n_times<10>(m1.col(0) +
m1.template triangularView<Lower>().solve(m1.col(0))),
2);
}
{
VERIFY(verify_eval_type<10>(m1, m1));
if (!NumTraits<Scalar>::IsComplex) {
VERIFY(verify_eval_type<3>(2 * m1, 2 * m1));
VERIFY(verify_eval_type<4>(2 * m1, m1));
} else {
VERIFY(verify_eval_type<2>(2 * m1, 2 * m1));
VERIFY(verify_eval_type<3>(2 * m1, m1));
}
VERIFY(verify_eval_type<2>(m1 + m1, m1 + m1));
VERIFY(verify_eval_type<3>(m1 + m1, m1));
VERIFY(verify_eval_type<1>(m1 * m1.transpose(), m2));
VERIFY(verify_eval_type<1>(m1 * (m1 + m1).transpose(), m2));
VERIFY(verify_eval_type<2>(m1 * m1.transpose(), m2));
VERIFY(verify_eval_type<1>(m1 + m1 * m1, m1));
VERIFY(
verify_eval_type<1>(m1.template triangularView<Lower>().solve(m1), m1));
VERIFY(verify_eval_type<1>(
m1 + m1.template triangularView<Lower>().solve(m1), m1));
}
}
EIGEN_DECLARE_TEST(nesting_ops) {
CALL_SUBTEST_1(run_nesting_ops_1(MatrixXf::Random(25, 25)));
CALL_SUBTEST_2(run_nesting_ops_1(MatrixXcd::Random(25, 25)));
CALL_SUBTEST_3(run_nesting_ops_1(Matrix4f::Random()));
CALL_SUBTEST_4(run_nesting_ops_1(Matrix2d::Random()));
Index s = internal::random<int>(1, EIGEN_TEST_MAX_SIZE);
CALL_SUBTEST_1(run_nesting_ops_2(MatrixXf(s, s)));
CALL_SUBTEST_2(run_nesting_ops_2(MatrixXcd(s, s)));
CALL_SUBTEST_3(run_nesting_ops_2(Matrix4f()));
CALL_SUBTEST_4(run_nesting_ops_2(Matrix2d()));
TEST_SET_BUT_UNUSED_VARIABLE(s)
}