aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/redux.cpp
blob: dee9d4bbd377e8686de01809199de6b067703e34 (plain)
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
122
123
124
125
126
127
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra. Eigen itself is part of the KDE project.
//
// Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 3 of the License, or (at your option) any later version.
//
// Alternatively, you can redistribute it and/or
// modify it under the terms of the GNU General Public License as
// published by the Free Software Foundation; either version 2 of
// the License, or (at your option) any later version.
//
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.

#include "main.h"

template<typename MatrixType> void matrixRedux(const MatrixType& m)
{
  typedef typename MatrixType::Scalar Scalar;

  int rows = m.rows();
  int cols = m.cols();

  MatrixType m1 = MatrixType::Random(rows, cols);

  VERIFY_IS_MUCH_SMALLER_THAN(MatrixType::Zero(rows, cols).sum(), Scalar(1));
  VERIFY_IS_APPROX(MatrixType::Ones(rows, cols).sum(), Scalar(float(rows*cols))); // the float() here to shut up excessive MSVC warning about int->complex conversion being lossy
  Scalar s(0), p(1), minc(ei_real(m1.coeff(0))), maxc(ei_real(m1.coeff(0)));
  for(int j = 0; j < cols; j++)
  for(int i = 0; i < rows; i++)
  {
    s += m1(i,j);
    p *= m1(i,j);
    minc = std::min(ei_real(minc), ei_real(m1(i,j)));
    maxc = std::max(ei_real(maxc), ei_real(m1(i,j)));
  }
  VERIFY_IS_APPROX(m1.sum(), s);
  VERIFY_IS_APPROX(m1.prod(), p);
  VERIFY_IS_APPROX(m1.real().minCoeff(), ei_real(minc));
  VERIFY_IS_APPROX(m1.real().maxCoeff(), ei_real(maxc));
}

template<typename VectorType> void vectorRedux(const VectorType& w)
{
  typedef typename VectorType::Scalar Scalar;
  typedef typename NumTraits<Scalar>::Real RealScalar;
  int size = w.size();

  VectorType v = VectorType::Random(size);
  for(int i = 1; i < size; i++)
  {
    Scalar s(0), p(1);
    RealScalar minc(ei_real(v.coeff(0))), maxc(ei_real(v.coeff(0)));
    for(int j = 0; j < i; j++)
    {
      s += v[j];
      p *= v[j];
      minc = std::min(minc, ei_real(v[j]));
      maxc = std::max(maxc, ei_real(v[j]));
    }
    VERIFY_IS_APPROX(s, v.start(i).sum());
    VERIFY_IS_APPROX(p, v.start(i).prod());
    VERIFY_IS_APPROX(minc, v.real().start(i).minCoeff());
    VERIFY_IS_APPROX(maxc, v.real().start(i).maxCoeff());
  }

  for(int i = 0; i < size-1; i++)
  {
    Scalar s(0), p(1);
    RealScalar minc(ei_real(v.coeff(i))), maxc(ei_real(v.coeff(i)));
    for(int j = i; j < size; j++)
    {
      s += v[j];
      p *= v[j];
      minc = std::min(minc, ei_real(v[j]));
      maxc = std::max(maxc, ei_real(v[j]));
    }
    VERIFY_IS_APPROX(s, v.end(size-i).sum());
    VERIFY_IS_APPROX(p, v.end(size-i).prod());
    VERIFY_IS_APPROX(minc, v.real().end(size-i).minCoeff());
    VERIFY_IS_APPROX(maxc, v.real().end(size-i).maxCoeff());
  }

  for(int i = 0; i < size/2; i++)
  {
    Scalar s(0), p(1);
    RealScalar minc(ei_real(v.coeff(i))), maxc(ei_real(v.coeff(i)));
    for(int j = i; j < size-i; j++)
    {
      s += v[j];
      p *= v[j];
      minc = std::min(minc, ei_real(v[j]));
      maxc = std::max(maxc, ei_real(v[j]));
    }
    VERIFY_IS_APPROX(s, v.segment(i, size-2*i).sum());
    VERIFY_IS_APPROX(p, v.segment(i, size-2*i).prod());
    VERIFY_IS_APPROX(minc, v.real().segment(i, size-2*i).minCoeff());
    VERIFY_IS_APPROX(maxc, v.real().segment(i, size-2*i).maxCoeff());
  }
}

void test_redux()
{
  for(int i = 0; i < g_repeat; i++) {
    CALL_SUBTEST( matrixRedux(Matrix<float, 1, 1>()) );
    CALL_SUBTEST( matrixRedux(Matrix2f()) );
    CALL_SUBTEST( matrixRedux(Matrix4d()) );
    CALL_SUBTEST( matrixRedux(MatrixXcf(3, 3)) );
    CALL_SUBTEST( matrixRedux(MatrixXd(8, 12)) );
    CALL_SUBTEST( matrixRedux(MatrixXi(8, 12)) );
  }
  for(int i = 0; i < g_repeat; i++) {
    CALL_SUBTEST( vectorRedux(VectorXf(5)) );
    CALL_SUBTEST( vectorRedux(VectorXd(10)) );
    CALL_SUBTEST( vectorRedux(VectorXf(33)) );
  }
}