blob: cd0f65f26d13bb35e57de8cf1976995c6a18a6a7 [file] [log] [blame]
Gael Guennebaud0be89a42009-03-05 10:25:22 +00001// This file is part of Eigen, a lightweight C++ template library
Benoit Jacob6347b1d2009-05-22 20:25:33 +02002// for linear algebra.
Gael Guennebaud0be89a42009-03-05 10:25:22 +00003//
4// Copyright (C) 2009 Gael Guennebaud <g.gael@free.fr>
5//
6// Eigen is free software; you can redistribute it and/or
7// modify it under the terms of the GNU Lesser General Public
8// License as published by the Free Software Foundation; either
9// version 3 of the License, or (at your option) any later version.
10//
11// Alternatively, you can redistribute it and/or
12// modify it under the terms of the GNU General Public License as
13// published by the Free Software Foundation; either version 2 of
14// the License, or (at your option) any later version.
15//
16// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
17// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
19// GNU General Public License for more details.
20//
21// You should have received a copy of the GNU Lesser General Public
22// License and a copy of the GNU General Public License along with
23// Eigen. If not, see <http://www.gnu.org/licenses/>.
24
25#include "main.h"
26#include <Eigen/Array>
27
28template<typename MatrixType> void replicate(const MatrixType& m)
29{
30 /* this test covers the following files:
31 Replicate.cpp
32 */
33
34 typedef typename MatrixType::Scalar Scalar;
35 typedef typename NumTraits<Scalar>::Real RealScalar;
36 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
37 typedef Matrix<Scalar, Dynamic, Dynamic> MatrixX;
38 typedef Matrix<Scalar, Dynamic, 1> VectorX;
39
40 int rows = m.rows();
41 int cols = m.cols();
42
43 MatrixType m1 = MatrixType::Random(rows, cols),
44 m2 = MatrixType::Random(rows, cols);
Gael Guennebaud14430942009-10-13 09:23:09 +020045
Gael Guennebaud0be89a42009-03-05 10:25:22 +000046 VectorType v1 = VectorType::Random(rows);
Gael Guennebaud14430942009-10-13 09:23:09 +020047
Gael Guennebaud0be89a42009-03-05 10:25:22 +000048 MatrixX x1, x2;
49 VectorX vx1;
50
51 int f1 = ei_random<int>(1,10),
52 f2 = ei_random<int>(1,10);
53
54 x1.resize(rows*f1,cols*f2);
55 for(int j=0; j<f2; j++)
56 for(int i=0; i<f1; i++)
57 x1.block(i*rows,j*cols,rows,cols) = m1;
58 VERIFY_IS_APPROX(x1, m1.replicate(f1,f2));
Gael Guennebaud14430942009-10-13 09:23:09 +020059
Gael Guennebaud0be89a42009-03-05 10:25:22 +000060 x2.resize(2*rows,3*cols);
61 x2 << m2, m2, m2,
62 m2, m2, m2;
63 VERIFY_IS_APPROX(x2, (m2.template replicate<2,3>()));
Gael Guennebaud14430942009-10-13 09:23:09 +020064
Gael Guennebaud0be89a42009-03-05 10:25:22 +000065 x2.resize(rows,f1);
66 for (int j=0; j<f1; ++j)
67 x2.col(j) = v1;
68 VERIFY_IS_APPROX(x2, v1.rowwise().replicate(f1));
Gael Guennebaud14430942009-10-13 09:23:09 +020069
Gael Guennebaud0be89a42009-03-05 10:25:22 +000070 vx1.resize(rows*f2);
71 for (int j=0; j<f2; ++j)
72 vx1.segment(j*rows,rows) = v1;
73 VERIFY_IS_APPROX(vx1, v1.colwise().replicate(f2));
74}
75
76void test_array_replicate()
77{
78 for(int i = 0; i < g_repeat; i++) {
79 CALL_SUBTEST( replicate(Matrix<float, 1, 1>()) );
80 CALL_SUBTEST( replicate(Vector2f()) );
81 CALL_SUBTEST( replicate(Vector3d()) );
82 CALL_SUBTEST( replicate(Vector4f()) );
83 CALL_SUBTEST( replicate(VectorXf(16)) );
84 CALL_SUBTEST( replicate(VectorXcd(10)) );
85 }
86}