blob: 1a917192b967f4bd6bb2e9d3da08c491ef72e14f [file] [log] [blame]
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +00001// This file is part of Eigen, a lightweight C++ template library
Benoit Jacob6347b1d2009-05-22 20:25:33 +02002// for linear algebra.
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +00003//
Gael Guennebaudf52d1192008-09-03 00:32:56 +00004// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
Benoit Jacob00f89a82008-11-24 13:40:43 +00005// Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +00006//
7// Eigen is free software; you can redistribute it and/or
8// modify it under the terms of the GNU Lesser General Public
9// License as published by the Free Software Foundation; either
10// version 3 of the License, or (at your option) any later version.
11//
12// Alternatively, you can redistribute it and/or
13// modify it under the terms of the GNU General Public License as
14// published by the Free Software Foundation; either version 2 of
15// the License, or (at your option) any later version.
16//
17// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
18// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
20// GNU General Public License for more details.
21//
22// You should have received a copy of the GNU Lesser General Public
23// License and a copy of the GNU General Public License along with
24// Eigen. If not, see <http://www.gnu.org/licenses/>.
25
Gael Guennebaude3fac692008-06-08 15:03:23 +000026// this hack is needed to make this file compiles with -pedantic (gcc)
Benoit Jacob9e3c7312009-01-22 00:09:34 +000027#ifdef __GNUC__
Gael Guennebaud60726f92008-06-04 10:15:48 +000028#define throw(X)
Benoit Jacob9e3c7312009-01-22 00:09:34 +000029#endif
Benoit Jacob1d52bd42009-01-08 15:20:21 +000030// discard stack allocation as that too bypasses malloc
Benoit Jacobd6712052008-12-31 00:25:31 +000031#define EIGEN_STACK_ALLOCATION_LIMIT 0
Benoit Jacob1d52bd42009-01-08 15:20:21 +000032// any heap allocation will raise an assert
33#define EIGEN_NO_MALLOC
Benoit Jacobd6712052008-12-31 00:25:31 +000034
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000035#include "main.h"
36
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000037template<typename MatrixType> void nomalloc(const MatrixType& m)
38{
39 /* this test check no dynamic memory allocation are issued with fixed-size matrices
40 */
41
42 typedef typename MatrixType::Scalar Scalar;
43 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
44
45 int rows = m.rows();
46 int cols = m.cols();
47
Gael Guennebaudc10f0692008-07-21 00:34:46 +000048 MatrixType m1 = MatrixType::Random(rows, cols),
49 m2 = MatrixType::Random(rows, cols),
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000050 m3(rows, cols),
Gael Guennebaudc10f0692008-07-21 00:34:46 +000051 mzero = MatrixType::Zero(rows, cols),
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000052 identity = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
Gael Guennebaudc10f0692008-07-21 00:34:46 +000053 ::Identity(rows, rows),
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000054 square = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
Gael Guennebaudc10f0692008-07-21 00:34:46 +000055 ::Random(rows, rows);
56 VectorType v1 = VectorType::Random(rows),
57 v2 = VectorType::Random(rows),
58 vzero = VectorType::Zero(rows);
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000059
60 Scalar s1 = ei_random<Scalar>();
61
62 int r = ei_random<int>(0, rows-1),
63 c = ei_random<int>(0, cols-1);
64
65 VERIFY_IS_APPROX((m1+m2)*s1, s1*m1+s1*m2);
66 VERIFY_IS_APPROX((m1+m2)(r,c), (m1(r,c))+(m2(r,c)));
Benoit Jacobf5791ee2008-07-08 00:49:10 +000067 VERIFY_IS_APPROX(m1.cwise() * m1.block(0,0,rows,cols), m1.cwise() * m1);
Gael Guennebaud56d00772009-08-06 12:20:02 +020068 if (MatrixType::RowsAtCompileTime<EIGEN_CACHEFRIENDLY_PRODUCT_THRESHOLD) {
69 // If the matrices are too large, we have better to use the optimized GEMM
70 // routines which allocates temporaries. However, on some platforms
71 // these temporaries are allocated on the stack using alloca.
72 VERIFY_IS_APPROX((m1*m1.transpose())*m2, m1*(m1.transpose()*m2));
73 }
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000074}
75
76void test_nomalloc()
77{
78 // check that our operator new is indeed called:
Gael Guennebaudc10f0692008-07-21 00:34:46 +000079 VERIFY_RAISES_ASSERT(MatrixXd dummy = MatrixXd::Random(3,3));
Benoit Jacob2840ac72009-10-28 18:19:29 -040080 CALL_SUBTEST(nomalloc(Matrix<float, 1, 1>()) );
81 CALL_SUBTEST(nomalloc(Matrix4d()) );
82 CALL_SUBTEST(nomalloc(Matrix<float,32,32>()) );
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000083}