blob: 8c4845d3cd16b223a1ac118f868f9f861eb671ee [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"
Gael Guennebaud0be89a42009-03-05 10:25:22 +000026
27template<typename MatrixType> void replicate(const MatrixType& m)
28{
29 /* this test covers the following files:
30 Replicate.cpp
31 */
32
33 typedef typename MatrixType::Scalar Scalar;
34 typedef typename NumTraits<Scalar>::Real RealScalar;
35 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
36 typedef Matrix<Scalar, Dynamic, Dynamic> MatrixX;
37 typedef Matrix<Scalar, Dynamic, 1> VectorX;
38
39 int rows = m.rows();
40 int cols = m.cols();
41
42 MatrixType m1 = MatrixType::Random(rows, cols),
43 m2 = MatrixType::Random(rows, cols);
Gael Guennebaud14430942009-10-13 09:23:09 +020044
Gael Guennebaud0be89a42009-03-05 10:25:22 +000045 VectorType v1 = VectorType::Random(rows);
Gael Guennebaud14430942009-10-13 09:23:09 +020046
Gael Guennebaud0be89a42009-03-05 10:25:22 +000047 MatrixX x1, x2;
48 VectorX vx1;
49
50 int f1 = ei_random<int>(1,10),
51 f2 = ei_random<int>(1,10);
52
53 x1.resize(rows*f1,cols*f2);
54 for(int j=0; j<f2; j++)
55 for(int i=0; i<f1; i++)
56 x1.block(i*rows,j*cols,rows,cols) = m1;
57 VERIFY_IS_APPROX(x1, m1.replicate(f1,f2));
Gael Guennebaud14430942009-10-13 09:23:09 +020058
Gael Guennebaud0be89a42009-03-05 10:25:22 +000059 x2.resize(2*rows,3*cols);
60 x2 << m2, m2, m2,
61 m2, m2, m2;
62 VERIFY_IS_APPROX(x2, (m2.template replicate<2,3>()));
Gael Guennebaud14430942009-10-13 09:23:09 +020063
Gael Guennebaud0be89a42009-03-05 10:25:22 +000064 x2.resize(rows,f1);
65 for (int j=0; j<f1; ++j)
66 x2.col(j) = v1;
67 VERIFY_IS_APPROX(x2, v1.rowwise().replicate(f1));
Gael Guennebaud14430942009-10-13 09:23:09 +020068
Gael Guennebaud0be89a42009-03-05 10:25:22 +000069 vx1.resize(rows*f2);
70 for (int j=0; j<f2; ++j)
71 vx1.segment(j*rows,rows) = v1;
72 VERIFY_IS_APPROX(vx1, v1.colwise().replicate(f2));
73}
74
75void test_array_replicate()
76{
77 for(int i = 0; i < g_repeat; i++) {
Benoit Jacob2840ac72009-10-28 18:19:29 -040078 CALL_SUBTEST_1( replicate(Matrix<float, 1, 1>()) );
79 CALL_SUBTEST_2( replicate(Vector2f()) );
80 CALL_SUBTEST_3( replicate(Vector3d()) );
81 CALL_SUBTEST_4( replicate(Vector4f()) );
82 CALL_SUBTEST_5( replicate(VectorXf(16)) );
83 CALL_SUBTEST_6( replicate(VectorXcd(10)) );
Gael Guennebaud0be89a42009-03-05 10:25:22 +000084 }
85}