2008-04-26 07:13:20 +08:00
|
|
|
// This file is part of Eigen, a lightweight C++ template library
|
2009-05-23 02:25:33 +08:00
|
|
|
// for linear algebra.
|
2008-04-26 07:13:20 +08:00
|
|
|
//
|
2008-11-24 21:40:43 +08:00
|
|
|
// Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
|
2010-06-25 05:21:58 +08:00
|
|
|
// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
|
2008-04-26 07:13:20 +08:00
|
|
|
//
|
2012-07-14 02:42:47 +08:00
|
|
|
// This Source Code Form is subject to the terms of the Mozilla
|
|
|
|
// Public License v. 2.0. If a copy of the MPL was not distributed
|
|
|
|
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
|
2008-04-26 07:13:20 +08:00
|
|
|
|
|
|
|
#include "main.h"
|
|
|
|
#include <Eigen/LU>
|
|
|
|
|
2008-08-04 12:45:59 +08:00
|
|
|
template<typename MatrixType> void determinant(const MatrixType& m)
|
2008-04-26 07:13:20 +08:00
|
|
|
{
|
|
|
|
/* this test covers the following files:
|
|
|
|
Determinant.h
|
|
|
|
*/
|
2010-06-21 17:36:00 +08:00
|
|
|
typedef typename MatrixType::Index Index;
|
|
|
|
Index size = m.rows();
|
2008-04-26 07:13:20 +08:00
|
|
|
|
2008-08-04 12:45:59 +08:00
|
|
|
MatrixType m1(size, size), m2(size, size);
|
|
|
|
m1.setRandom();
|
|
|
|
m2.setRandom();
|
2008-04-26 07:13:20 +08:00
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
2010-10-25 22:15:22 +08:00
|
|
|
Scalar x = internal::random<Scalar>();
|
2009-01-04 23:26:32 +08:00
|
|
|
VERIFY_IS_APPROX(MatrixType::Identity(size, size).determinant(), Scalar(1));
|
2009-08-06 22:54:55 +08:00
|
|
|
VERIFY_IS_APPROX((m1*m2).eval().determinant(), m1.determinant() * m2.determinant());
|
2008-08-04 12:45:59 +08:00
|
|
|
if(size==1) return;
|
2010-10-25 22:15:22 +08:00
|
|
|
Index i = internal::random<Index>(0, size-1);
|
2010-06-21 17:36:00 +08:00
|
|
|
Index j;
|
2008-08-04 12:45:59 +08:00
|
|
|
do {
|
2010-10-25 22:15:22 +08:00
|
|
|
j = internal::random<Index>(0, size-1);
|
2008-08-04 12:45:59 +08:00
|
|
|
} while(j==i);
|
|
|
|
m2 = m1;
|
|
|
|
m2.row(i).swap(m2.row(j));
|
2009-01-04 23:26:32 +08:00
|
|
|
VERIFY_IS_APPROX(m2.determinant(), -m1.determinant());
|
2008-08-04 12:45:59 +08:00
|
|
|
m2 = m1;
|
|
|
|
m2.col(i).swap(m2.col(j));
|
2009-01-04 23:26:32 +08:00
|
|
|
VERIFY_IS_APPROX(m2.determinant(), -m1.determinant());
|
|
|
|
VERIFY_IS_APPROX(m2.determinant(), m2.transpose().determinant());
|
2013-06-11 05:40:56 +08:00
|
|
|
VERIFY_IS_APPROX(numext::conj(m2.determinant()), m2.adjoint().determinant());
|
2008-08-04 12:45:59 +08:00
|
|
|
m2 = m1;
|
|
|
|
m2.row(i) += x*m2.row(j);
|
2009-01-04 23:26:32 +08:00
|
|
|
VERIFY_IS_APPROX(m2.determinant(), m1.determinant());
|
2008-08-04 12:45:59 +08:00
|
|
|
m2 = m1;
|
|
|
|
m2.row(i) *= x;
|
2009-01-04 23:26:32 +08:00
|
|
|
VERIFY_IS_APPROX(m2.determinant(), m1.determinant() * x);
|
2010-07-19 16:45:06 +08:00
|
|
|
|
|
|
|
// check empty matrix
|
|
|
|
VERIFY_IS_APPROX(m2.block(0,0,0,0).determinant(), Scalar(1));
|
2008-04-26 07:13:20 +08:00
|
|
|
}
|
|
|
|
|
2008-05-22 20:18:55 +08:00
|
|
|
void test_determinant()
|
2008-04-26 07:13:20 +08:00
|
|
|
{
|
2008-05-22 20:18:55 +08:00
|
|
|
for(int i = 0; i < g_repeat; i++) {
|
2013-06-24 01:11:32 +08:00
|
|
|
int s = 0;
|
2009-10-29 06:19:29 +08:00
|
|
|
CALL_SUBTEST_1( determinant(Matrix<float, 1, 1>()) );
|
|
|
|
CALL_SUBTEST_2( determinant(Matrix<double, 2, 2>()) );
|
|
|
|
CALL_SUBTEST_3( determinant(Matrix<double, 3, 3>()) );
|
|
|
|
CALL_SUBTEST_4( determinant(Matrix<double, 4, 4>()) );
|
|
|
|
CALL_SUBTEST_5( determinant(Matrix<std::complex<double>, 10, 10>()) );
|
2011-07-12 20:41:00 +08:00
|
|
|
s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/4);
|
|
|
|
CALL_SUBTEST_6( determinant(MatrixXd(s, s)) );
|
2013-06-25 17:42:04 +08:00
|
|
|
TEST_SET_BUT_UNUSED_VARIABLE(s)
|
2008-04-26 07:13:20 +08:00
|
|
|
}
|
|
|
|
}
|