aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/array_replicate.cpp
blob: d1608915f5a69c2b0f203192b141c7c422979815 (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
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2009 Gael Guennebaud <g.gael@free.fr>
//
// 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"
#include <Eigen/Array>

template<typename MatrixType> void replicate(const MatrixType& m)
{
  /* this test covers the following files:
     Replicate.cpp
  */

  typedef typename MatrixType::Scalar Scalar;
  typedef typename NumTraits<Scalar>::Real RealScalar;
  typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  typedef Matrix<Scalar, Dynamic, Dynamic> MatrixX;
  typedef Matrix<Scalar, Dynamic, 1> VectorX;

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

  MatrixType m1 = MatrixType::Random(rows, cols),
             m2 = MatrixType::Random(rows, cols);
  
  VectorType v1 = VectorType::Random(rows);
  
  MatrixX x1, x2;
  VectorX vx1;

  int  f1 = ei_random<int>(1,10),
       f2 = ei_random<int>(1,10);

  x1.resize(rows*f1,cols*f2);
  for(int j=0; j<f2; j++)
  for(int i=0; i<f1; i++)
    x1.block(i*rows,j*cols,rows,cols) = m1;
  VERIFY_IS_APPROX(x1, m1.replicate(f1,f2));
  
  x2.resize(2*rows,3*cols);
  x2 << m2, m2, m2,
        m2, m2, m2;
  VERIFY_IS_APPROX(x2, (m2.template replicate<2,3>()));
  
  x2.resize(rows,f1);
  for (int j=0; j<f1; ++j)
    x2.col(j) = v1;
  VERIFY_IS_APPROX(x2, v1.rowwise().replicate(f1));
  
  vx1.resize(rows*f2);
  for (int j=0; j<f2; ++j)
    vx1.segment(j*rows,rows) = v1;
  VERIFY_IS_APPROX(vx1, v1.colwise().replicate(f2));
}

void test_array_replicate()
{
  for(int i = 0; i < g_repeat; i++) {
    CALL_SUBTEST( replicate(Matrix<float, 1, 1>()) );
    CALL_SUBTEST( replicate(Vector2f()) );
    CALL_SUBTEST( replicate(Vector3d()) );
    CALL_SUBTEST( replicate(Vector4f()) );
    CALL_SUBTEST( replicate(VectorXf(16)) );
    CALL_SUBTEST( replicate(VectorXcd(10)) );
  }
}