blob: 94da7425b87a1327fc2df2bb3f459d2444636228 [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//
Gael Guennebaud28e64b02010-06-24 23:21:58 +02004// Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
Gael Guennebaud0be89a42009-03-05 10:25:22 +00005//
Benoit Jacob69124cf2012-07-13 14:42:47 -04006// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
Gael Guennebaud0be89a42009-03-05 10:25:22 +00009
10#include "main.h"
Gael Guennebaud0be89a42009-03-05 10:25:22 +000011
12template<typename MatrixType> void replicate(const MatrixType& m)
13{
14 /* this test covers the following files:
15 Replicate.cpp
16 */
Hauke Heibelf1679c72010-06-20 17:37:56 +020017 typedef typename MatrixType::Index Index;
Gael Guennebaud0be89a42009-03-05 10:25:22 +000018 typedef typename MatrixType::Scalar Scalar;
19 typedef typename NumTraits<Scalar>::Real RealScalar;
20 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
21 typedef Matrix<Scalar, Dynamic, Dynamic> MatrixX;
22 typedef Matrix<Scalar, Dynamic, 1> VectorX;
23
Hauke Heibelf1679c72010-06-20 17:37:56 +020024 Index rows = m.rows();
25 Index cols = m.cols();
Gael Guennebaud0be89a42009-03-05 10:25:22 +000026
27 MatrixType m1 = MatrixType::Random(rows, cols),
28 m2 = MatrixType::Random(rows, cols);
Gael Guennebaud14430942009-10-13 09:23:09 +020029
Gael Guennebaud0be89a42009-03-05 10:25:22 +000030 VectorType v1 = VectorType::Random(rows);
Gael Guennebaud14430942009-10-13 09:23:09 +020031
Gael Guennebaud0be89a42009-03-05 10:25:22 +000032 MatrixX x1, x2;
33 VectorX vx1;
34
Benoit Jacob47160402010-10-25 10:15:22 -040035 int f1 = internal::random<int>(1,10),
36 f2 = internal::random<int>(1,10);
Gael Guennebaud0be89a42009-03-05 10:25:22 +000037
38 x1.resize(rows*f1,cols*f2);
39 for(int j=0; j<f2; j++)
40 for(int i=0; i<f1; i++)
41 x1.block(i*rows,j*cols,rows,cols) = m1;
42 VERIFY_IS_APPROX(x1, m1.replicate(f1,f2));
Gael Guennebaud14430942009-10-13 09:23:09 +020043
Gael Guennebaud0be89a42009-03-05 10:25:22 +000044 x2.resize(2*rows,3*cols);
45 x2 << m2, m2, m2,
46 m2, m2, m2;
47 VERIFY_IS_APPROX(x2, (m2.template replicate<2,3>()));
Gael Guennebaud14430942009-10-13 09:23:09 +020048
Gael Guennebaud0be89a42009-03-05 10:25:22 +000049 x2.resize(rows,f1);
50 for (int j=0; j<f1; ++j)
51 x2.col(j) = v1;
52 VERIFY_IS_APPROX(x2, v1.rowwise().replicate(f1));
Gael Guennebaud14430942009-10-13 09:23:09 +020053
Gael Guennebaud0be89a42009-03-05 10:25:22 +000054 vx1.resize(rows*f2);
55 for (int j=0; j<f2; ++j)
56 vx1.segment(j*rows,rows) = v1;
57 VERIFY_IS_APPROX(vx1, v1.colwise().replicate(f2));
58}
59
60void test_array_replicate()
61{
62 for(int i = 0; i < g_repeat; i++) {
Benoit Jacob2840ac72009-10-28 18:19:29 -040063 CALL_SUBTEST_1( replicate(Matrix<float, 1, 1>()) );
64 CALL_SUBTEST_2( replicate(Vector2f()) );
65 CALL_SUBTEST_3( replicate(Vector3d()) );
66 CALL_SUBTEST_4( replicate(Vector4f()) );
67 CALL_SUBTEST_5( replicate(VectorXf(16)) );
68 CALL_SUBTEST_6( replicate(VectorXcd(10)) );
Gael Guennebaud0be89a42009-03-05 10:25:22 +000069 }
70}