mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-04-22 09:39:34 +08:00
* bug fixes in: Dot, generalized eigen problem, singular matrix detetection in Cholesky
* fix all numerical instabilies in the unit tests, now all tests can be run 2000 times with almost zero failures.
This commit is contained in:
parent
312013a089
commit
2120fed849
@ -93,17 +93,18 @@ void Cholesky<MatrixType>::compute(const MatrixType& a)
|
||||
assert(a.rows()==a.cols());
|
||||
const int size = a.rows();
|
||||
m_matrix.resize(size, size);
|
||||
const RealScalar eps = ei_sqrt(precision<Scalar>());
|
||||
|
||||
RealScalar x;
|
||||
x = ei_real(a.coeff(0,0));
|
||||
m_isPositiveDefinite = x > precision<Scalar>() && ei_isMuchSmallerThan(ei_imag(a.coeff(0,0)), RealScalar(1));
|
||||
m_isPositiveDefinite = x > eps && ei_isMuchSmallerThan(ei_imag(a.coeff(0,0)), RealScalar(1));
|
||||
m_matrix.coeffRef(0,0) = ei_sqrt(x);
|
||||
m_matrix.col(0).end(size-1) = a.row(0).end(size-1).adjoint() / ei_real(m_matrix.coeff(0,0));
|
||||
for (int j = 1; j < size; ++j)
|
||||
{
|
||||
Scalar tmp = ei_real(a.coeff(j,j)) - m_matrix.row(j).start(j).norm2();
|
||||
x = ei_real(tmp);
|
||||
if (x < precision<Scalar>() || (!ei_isMuchSmallerThan(ei_imag(tmp), RealScalar(1))))
|
||||
if (x < eps || (!ei_isMuchSmallerThan(ei_imag(tmp), RealScalar(1))))
|
||||
{
|
||||
m_isPositiveDefinite = false;
|
||||
return;
|
||||
|
@ -94,6 +94,7 @@ void CholeskyWithoutSquareRoot<MatrixType>::compute(const MatrixType& a)
|
||||
const int size = a.rows();
|
||||
m_matrix.resize(size, size);
|
||||
m_isPositiveDefinite = true;
|
||||
const RealScalar eps = ei_sqrt(precision<Scalar>());
|
||||
|
||||
// Let's preallocate a temporay vector to evaluate the matrix-vector product into it.
|
||||
// Unlike the standard Cholesky decomposition, here we cannot evaluate it to the destination
|
||||
@ -111,7 +112,7 @@ void CholeskyWithoutSquareRoot<MatrixType>::compute(const MatrixType& a)
|
||||
RealScalar tmp = ei_real(a.coeff(j,j) - (m_matrix.row(j).start(j) * m_matrix.col(j).start(j).conjugate()).coeff(0,0));
|
||||
m_matrix.coeffRef(j,j) = tmp;
|
||||
|
||||
if (ei_isMuchSmallerThan(tmp,RealScalar(1)))
|
||||
if (tmp < eps)
|
||||
{
|
||||
m_isPositiveDefinite = false;
|
||||
return;
|
||||
|
@ -229,9 +229,9 @@ struct ei_dot_impl<Derived1, Derived2, LinearVectorization, CompleteUnrolling>
|
||||
};
|
||||
static Scalar run(const Derived1& v1, const Derived2& v2)
|
||||
{
|
||||
Scalar res = ei_predux(ei_dot_vec_unroller<Derived1, Derived2, 0, VectorizationSize>::run(v1, v2));
|
||||
Scalar res = ei_predux(ei_dot_vec_unroller<Derived1, Derived2, 0, VectorizationSize>::run(v1, v2));
|
||||
if (VectorizationSize != Size)
|
||||
res += ei_dot_novec_unroller<Derived1, Derived2, VectorizationSize, Size>::run(v1, v2);
|
||||
res += ei_dot_novec_unroller<Derived1, Derived2, VectorizationSize, Size-VectorizationSize>::run(v1, v2);
|
||||
return res;
|
||||
}
|
||||
};
|
||||
|
@ -131,7 +131,7 @@ template<typename Scalar>
|
||||
AngleAxis<Scalar>& AngleAxis<Scalar>::operator=(const QuaternionType& q)
|
||||
{
|
||||
Scalar n2 = q.vec().norm2();
|
||||
if (ei_isMuchSmallerThan(n2,Scalar(1)))
|
||||
if (n2 < precision<Scalar>()*precision<Scalar>())
|
||||
{
|
||||
m_angle = 0;
|
||||
m_axis << 1, 0, 0;
|
||||
|
@ -226,21 +226,32 @@ compute(const MatrixType& matA, const MatrixType& matB, bool computeEigenvectors
|
||||
{
|
||||
ei_assert(matA.cols()==matA.rows() && matB.rows()==matA.rows() && matB.cols()==matB.rows());
|
||||
|
||||
// Compute the cholesky decomposition of matB = U'U
|
||||
// Compute the cholesky decomposition of matB = L L'
|
||||
Cholesky<MatrixType> cholB(matB);
|
||||
|
||||
// compute C = inv(U') A inv(U)
|
||||
MatrixType matC = cholB.matrixL().solveTriangular(matA);
|
||||
// FIXME since we currently do not support A * inv(U),
|
||||
// let's do (inv(U') A')' :
|
||||
matC = (cholB.matrixL().solveTriangular(matC.adjoint())).adjoint();
|
||||
// compute C = inv(L) A inv(L')
|
||||
MatrixType matC = matA;
|
||||
cholB.matrixL().solveTriangularInPlace(matC);
|
||||
// FIXME since we currently do not support A * inv(L'), let's do (inv(L) A')' :
|
||||
matC = matC.adjoint().eval();
|
||||
cholB.matrixL().template marked<Lower>().solveTriangularInPlace(matC);
|
||||
matC = matC.adjoint().eval();
|
||||
// this version works too:
|
||||
// matC = matC.transpose();
|
||||
// cholB.matrixL().conjugate().template marked<Lower>().solveTriangularInPlace(matC);
|
||||
// matC = matC.transpose();
|
||||
// FIXME: this should work: (currently it only does for small matrices)
|
||||
// Transpose<MatrixType> trMatC(matC);
|
||||
// cholB.matrixL().conjugate().eval().template marked<Lower>().solveTriangularInPlace(trMatC);
|
||||
|
||||
compute(matC, computeEigenvectors);
|
||||
|
||||
if (computeEigenvectors)
|
||||
{
|
||||
// transform back the eigen vectors: evecs = inv(U) * evecs
|
||||
m_eivec = cholB.matrixL().adjoint().template marked<Upper>().solveTriangular(m_eivec);
|
||||
cholB.matrixL().adjoint().template marked<Upper>().solveTriangularInPlace(m_eivec);
|
||||
for (int i=0; i<m_eivec.cols(); ++i)
|
||||
m_eivec.col(i) = m_eivec.col(i).normalized();
|
||||
}
|
||||
}
|
||||
|
||||
|
159
cmake/FindGSL.cmake
Normal file
159
cmake/FindGSL.cmake
Normal file
@ -0,0 +1,159 @@
|
||||
# Try to find gnu scientific library GSL
|
||||
# See
|
||||
# http://www.gnu.org/software/gsl/ and
|
||||
# http://gnuwin32.sourceforge.net/packages/gsl.htm
|
||||
#
|
||||
# Once run this will define:
|
||||
#
|
||||
# GSL_FOUND = system has GSL lib
|
||||
#
|
||||
# GSL_LIBRARIES = full path to the libraries
|
||||
# on Unix/Linux with additional linker flags from "gsl-config --libs"
|
||||
#
|
||||
# CMAKE_GSL_CXX_FLAGS = Unix compiler flags for GSL, essentially "`gsl-config --cxxflags`"
|
||||
#
|
||||
# GSL_INCLUDE_DIR = where to find headers
|
||||
#
|
||||
# GSL_LINK_DIRECTORIES = link directories, useful for rpath on Unix
|
||||
# GSL_EXE_LINKER_FLAGS = rpath on Unix
|
||||
#
|
||||
# Felix Woelk 07/2004
|
||||
# Jan Woetzel
|
||||
#
|
||||
# www.mip.informatik.uni-kiel.de
|
||||
# --------------------------------
|
||||
|
||||
IF(WIN32)
|
||||
# JW tested with gsl-1.8, Windows XP, MSVS 7.1
|
||||
SET(GSL_POSSIBLE_ROOT_DIRS
|
||||
${GSL_ROOT_DIR}
|
||||
$ENV{GSL_ROOT_DIR}
|
||||
${GSL_DIR}
|
||||
${GSL_HOME}
|
||||
$ENV{GSL_DIR}
|
||||
$ENV{GSL_HOME}
|
||||
$ENV{EXTRA}
|
||||
"C:/Program Files/GnuWin32"
|
||||
)
|
||||
FIND_PATH(GSL_INCLUDE_DIR
|
||||
NAMES gsl/gsl_cdf.h gsl/gsl_randist.h
|
||||
PATHS ${GSL_POSSIBLE_ROOT_DIRS}
|
||||
PATH_SUFFIXES include
|
||||
DOC "GSL header include dir"
|
||||
)
|
||||
|
||||
FIND_LIBRARY(GSL_GSL_LIBRARY
|
||||
NAMES libgsl.dll.a gsl libgsl
|
||||
PATHS ${GSL_POSSIBLE_ROOT_DIRS}
|
||||
PATH_SUFFIXES lib
|
||||
DOC "GSL library" )
|
||||
|
||||
if(NOT GSL_GSL_LIBRARY)
|
||||
FIND_FILE(GSL_GSL_LIBRARY
|
||||
NAMES libgsl.dll.a
|
||||
PATHS ${GSL_POSSIBLE_ROOT_DIRS}
|
||||
PATH_SUFFIXES lib
|
||||
DOC "GSL library")
|
||||
endif(NOT GSL_GSL_LIBRARY)
|
||||
|
||||
FIND_LIBRARY(GSL_GSLCBLAS_LIBRARY
|
||||
NAMES libgslcblas.dll.a gslcblas libgslcblas
|
||||
PATHS ${GSL_POSSIBLE_ROOT_DIRS}
|
||||
PATH_SUFFIXES lib
|
||||
DOC "GSL cblas library dir" )
|
||||
|
||||
if(NOT GSL_GSLCBLAS_LIBRARY)
|
||||
FIND_FILE(GSL_GSLCBLAS_LIBRARY
|
||||
NAMES libgslcblas.dll.a
|
||||
PATHS ${GSL_POSSIBLE_ROOT_DIRS}
|
||||
PATH_SUFFIXES lib
|
||||
DOC "GSL library")
|
||||
endif(NOT GSL_GSLCBLAS_LIBRARY)
|
||||
|
||||
SET(GSL_LIBRARIES ${GSL_GSL_LIBRARY})
|
||||
|
||||
#MESSAGE("DBG\n"
|
||||
# "GSL_GSL_LIBRARY=${GSL_GSL_LIBRARY}\n"
|
||||
# "GSL_GSLCBLAS_LIBRARY=${GSL_GSLCBLAS_LIBRARY}\n"
|
||||
# "GSL_LIBRARIES=${GSL_LIBRARIES}")
|
||||
|
||||
|
||||
ELSE(WIN32)
|
||||
|
||||
IF(UNIX)
|
||||
SET(GSL_CONFIG_PREFER_PATH
|
||||
"$ENV{GSL_DIR}/bin"
|
||||
"$ENV{GSL_DIR}"
|
||||
"$ENV{GSL_HOME}/bin"
|
||||
"$ENV{GSL_HOME}"
|
||||
CACHE STRING "preferred path to GSL (gsl-config)")
|
||||
FIND_PROGRAM(GSL_CONFIG gsl-config
|
||||
${GSL_CONFIG_PREFER_PATH}
|
||||
/usr/bin/
|
||||
)
|
||||
# MESSAGE("DBG GSL_CONFIG ${GSL_CONFIG}")
|
||||
|
||||
IF (GSL_CONFIG)
|
||||
# set CXXFLAGS to be fed into CXX_FLAGS by the user:
|
||||
SET(GSL_CXX_FLAGS "`${GSL_CONFIG} --cflags`")
|
||||
|
||||
# set INCLUDE_DIRS to prefix+include
|
||||
EXEC_PROGRAM(${GSL_CONFIG}
|
||||
ARGS --prefix
|
||||
OUTPUT_VARIABLE GSL_PREFIX)
|
||||
SET(GSL_INCLUDE_DIR ${GSL_PREFIX}/include CACHE STRING INTERNAL)
|
||||
|
||||
# set link libraries and link flags
|
||||
#SET(GSL_LIBRARIES "`${GSL_CONFIG} --libs`")
|
||||
EXEC_PROGRAM(${GSL_CONFIG}
|
||||
ARGS --libs
|
||||
OUTPUT_VARIABLE GSL_LIBRARIES )
|
||||
|
||||
# extract link dirs for rpath
|
||||
EXEC_PROGRAM(${GSL_CONFIG}
|
||||
ARGS --libs
|
||||
OUTPUT_VARIABLE GSL_CONFIG_LIBS )
|
||||
|
||||
# split off the link dirs (for rpath)
|
||||
# use regular expression to match wildcard equivalent "-L*<endchar>"
|
||||
# with <endchar> is a space or a semicolon
|
||||
STRING(REGEX MATCHALL "[-][L]([^ ;])+"
|
||||
GSL_LINK_DIRECTORIES_WITH_PREFIX
|
||||
"${GSL_CONFIG_LIBS}" )
|
||||
# MESSAGE("DBG GSL_LINK_DIRECTORIES_WITH_PREFIX=${GSL_LINK_DIRECTORIES_WITH_PREFIX}")
|
||||
|
||||
# remove prefix -L because we need the pure directory for LINK_DIRECTORIES
|
||||
|
||||
IF (GSL_LINK_DIRECTORIES_WITH_PREFIX)
|
||||
STRING(REGEX REPLACE "[-][L]" "" GSL_LINK_DIRECTORIES ${GSL_LINK_DIRECTORIES_WITH_PREFIX} )
|
||||
ENDIF (GSL_LINK_DIRECTORIES_WITH_PREFIX)
|
||||
SET(GSL_EXE_LINKER_FLAGS "-Wl,-rpath,${GSL_LINK_DIRECTORIES}" CACHE STRING INTERNAL)
|
||||
# MESSAGE("DBG GSL_LINK_DIRECTORIES=${GSL_LINK_DIRECTORIES}")
|
||||
# MESSAGE("DBG GSL_EXE_LINKER_FLAGS=${GSL_EXE_LINKER_FLAGS}")
|
||||
|
||||
# ADD_DEFINITIONS("-DHAVE_GSL")
|
||||
# SET(GSL_DEFINITIONS "-DHAVE_GSL")
|
||||
MARK_AS_ADVANCED(
|
||||
GSL_CXX_FLAGS
|
||||
GSL_INCLUDE_DIR
|
||||
GSL_LIBRARIES
|
||||
GSL_LINK_DIRECTORIES
|
||||
GSL_DEFINITIONS
|
||||
)
|
||||
MESSAGE(STATUS "Using GSL from ${GSL_PREFIX}")
|
||||
|
||||
ELSE(GSL_CONFIG)
|
||||
MESSAGE("FindGSL.cmake: gsl-config not found. Please set it manually. GSL_CONFIG=${GSL_CONFIG}")
|
||||
ENDIF(GSL_CONFIG)
|
||||
|
||||
ENDIF(UNIX)
|
||||
ENDIF(WIN32)
|
||||
|
||||
|
||||
IF(GSL_LIBRARIES)
|
||||
IF(GSL_INCLUDE_DIR OR GSL_CXX_FLAGS)
|
||||
|
||||
SET(GSL_FOUND 1)
|
||||
|
||||
ENDIF(GSL_INCLUDE_DIR OR GSL_CXX_FLAGS)
|
||||
ENDIF(GSL_LIBRARIES)
|
@ -1,5 +1,9 @@
|
||||
IF(BUILD_TESTS)
|
||||
|
||||
find_package(GSL)
|
||||
if(GSL_FOUND)
|
||||
add_definitions("-DHAS_GSL")
|
||||
endif(GSL_FOUND)
|
||||
|
||||
IF(CMAKE_COMPILER_IS_GNUCXX)
|
||||
IF(CMAKE_SYSTEM_NAME MATCHES Linux)
|
||||
@ -69,6 +73,10 @@ MACRO(EI_ADD_TEST testname)
|
||||
target_link_libraries(${targetname} Eigen2)
|
||||
ENDIF(TEST_LIB)
|
||||
|
||||
if(GSL_FOUND)
|
||||
target_link_libraries(${targetname} ${GSL_LIBRARIES})
|
||||
endif(GSL_FOUND)
|
||||
|
||||
IF(WIN32)
|
||||
ADD_TEST(${testname} "${targetname}")
|
||||
ELSE(WIN32)
|
||||
|
@ -31,25 +31,29 @@ template<typename MatrixType> void adjoint(const MatrixType& m)
|
||||
*/
|
||||
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
|
||||
int rows = m.rows();
|
||||
int cols = m.cols();
|
||||
|
||||
MatrixType m1 = MatrixType::Random(rows, cols),
|
||||
m2 = MatrixType::Random(rows, cols),
|
||||
RealScalar largerEps = test_precision<RealScalar>();
|
||||
if (ei_is_same_type<RealScalar,float>::ret)
|
||||
largerEps = 1e-3f;
|
||||
|
||||
MatrixType m1 = test_random_matrix<MatrixType>(rows, cols),
|
||||
m2 = test_random_matrix<MatrixType>(rows, cols),
|
||||
m3(rows, cols),
|
||||
mzero = MatrixType::Zero(rows, cols),
|
||||
identity = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
::Identity(rows, rows),
|
||||
square = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
::Random(rows, rows);
|
||||
VectorType v1 = VectorType::Random(rows),
|
||||
v2 = VectorType::Random(rows),
|
||||
v3 = VectorType::Random(rows),
|
||||
identity = SquareMatrixType::Identity(rows, rows),
|
||||
square = test_random_matrix<SquareMatrixType>(rows, rows);
|
||||
VectorType v1 = test_random_matrix<VectorType>(rows),
|
||||
v2 = test_random_matrix<VectorType>(rows),
|
||||
v3 = test_random_matrix<VectorType>(rows),
|
||||
vzero = VectorType::Zero(rows);
|
||||
|
||||
Scalar s1 = ei_random<Scalar>(),
|
||||
s2 = ei_random<Scalar>();
|
||||
Scalar s1 = test_random<Scalar>(),
|
||||
s2 = test_random<Scalar>();
|
||||
|
||||
// check basic compatibility of adjoint, transpose, conjugate
|
||||
VERIFY_IS_APPROX(m1.transpose().conjugate().adjoint(), m1);
|
||||
@ -61,19 +65,18 @@ template<typename MatrixType> void adjoint(const MatrixType& m)
|
||||
|
||||
// check basic properties of dot, norm, norm2
|
||||
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||
VERIFY_IS_APPROX((s1 * v1 + s2 * v2).dot(v3), s1 * v1.dot(v3) + s2 * v2.dot(v3));
|
||||
VERIFY_IS_APPROX(v3.dot(s1 * v1 + s2 * v2), ei_conj(s1)*v3.dot(v1)+ei_conj(s2)*v3.dot(v2));
|
||||
VERIFY_IS_APPROX(ei_conj(v1.dot(v2)), v2.dot(v1));
|
||||
VERIFY_IS_APPROX(ei_abs(v1.dot(v1)), v1.norm2());
|
||||
VERIFY(ei_isApprox((s1 * v1 + s2 * v2).dot(v3), s1 * v1.dot(v3) + s2 * v2.dot(v3), largerEps));
|
||||
VERIFY(ei_isApprox(v3.dot(s1 * v1 + s2 * v2), ei_conj(s1)*v3.dot(v1)+ei_conj(s2)*v3.dot(v2), largerEps));
|
||||
VERIFY_IS_APPROX(ei_conj(v1.dot(v2)), v2.dot(v1));
|
||||
VERIFY_IS_APPROX(ei_abs(v1.dot(v1)), v1.norm2());
|
||||
if(NumTraits<Scalar>::HasFloatingPoint)
|
||||
VERIFY_IS_APPROX(v1.norm2(), v1.norm() * v1.norm());
|
||||
VERIFY_IS_MUCH_SMALLER_THAN(ei_abs(vzero.dot(v1)), static_cast<RealScalar>(1));
|
||||
VERIFY_IS_APPROX(v1.norm2(), v1.norm() * v1.norm());
|
||||
VERIFY_IS_MUCH_SMALLER_THAN(ei_abs(vzero.dot(v1)), static_cast<RealScalar>(1));
|
||||
if(NumTraits<Scalar>::HasFloatingPoint)
|
||||
VERIFY_IS_MUCH_SMALLER_THAN(vzero.norm(), static_cast<RealScalar>(1));
|
||||
VERIFY_IS_MUCH_SMALLER_THAN(vzero.norm(), static_cast<RealScalar>(1));
|
||||
|
||||
// check compatibility of dot and adjoint
|
||||
// FIXME this line failed with MSVC and complex<double> in the ei_aligned_free()
|
||||
VERIFY_IS_APPROX(v1.dot(square * v2), (square.adjoint() * v1).dot(v2));
|
||||
VERIFY(ei_isApprox(v1.dot(square * v2), (square.adjoint() * v1).dot(v2), largerEps));
|
||||
|
||||
// like in testBasicStuff, test operator() to check const-qualification
|
||||
int r = ei_random<int>(0, rows-1),
|
||||
@ -93,10 +96,11 @@ void test_adjoint()
|
||||
{
|
||||
for(int i = 0; i < g_repeat; i++) {
|
||||
CALL_SUBTEST( adjoint(Matrix<float, 1, 1>()) );
|
||||
CALL_SUBTEST( adjoint(Matrix4d()) );
|
||||
CALL_SUBTEST( adjoint(MatrixXcf(3, 3)) );
|
||||
CALL_SUBTEST( adjoint(Matrix3d()) );
|
||||
CALL_SUBTEST( adjoint(Matrix4f()) );
|
||||
CALL_SUBTEST( adjoint(MatrixXcf(4, 4)) );
|
||||
CALL_SUBTEST( adjoint(MatrixXi(8, 12)) );
|
||||
CALL_SUBTEST( adjoint(MatrixXcd(20, 20)) );
|
||||
CALL_SUBTEST( adjoint(MatrixXf(21, 21)) );
|
||||
}
|
||||
// test a large matrix only once
|
||||
CALL_SUBTEST( adjoint(Matrix<float, 100, 100>()) );
|
||||
|
@ -32,17 +32,18 @@ template<typename MatrixType> void scalarAdd(const MatrixType& m)
|
||||
*/
|
||||
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
||||
|
||||
int rows = m.rows();
|
||||
int cols = m.cols();
|
||||
|
||||
MatrixType m1 = MatrixType::Random(rows, cols),
|
||||
m2 = MatrixType::Random(rows, cols),
|
||||
MatrixType m1 = test_random_matrix<MatrixType>(rows, cols),
|
||||
m2 = test_random_matrix<MatrixType>(rows, cols),
|
||||
m3(rows, cols);
|
||||
|
||||
Scalar s1 = ei_random<Scalar>(),
|
||||
s2 = ei_random<Scalar>();
|
||||
Scalar s1 = test_random<Scalar>(),
|
||||
s2 = test_random<Scalar>();
|
||||
|
||||
VERIFY_IS_APPROX(m1.cwise() + s1, s1 + m1.cwise());
|
||||
VERIFY_IS_APPROX(m1.cwise() + s1, MatrixType::Constant(rows,cols,s1) + m1);
|
||||
@ -56,7 +57,8 @@ template<typename MatrixType> void scalarAdd(const MatrixType& m)
|
||||
|
||||
VERIFY_IS_APPROX(m1.colwise().sum().sum(), m1.sum());
|
||||
VERIFY_IS_APPROX(m1.rowwise().sum().sum(), m1.sum());
|
||||
VERIFY_IS_NOT_APPROX((m1.rowwise().sum()*2).sum(), m1.sum());
|
||||
if (!ei_isApprox(m1.sum(), (m1+m2).sum()))
|
||||
VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
|
||||
VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(ei_scalar_sum_op<Scalar>()));
|
||||
}
|
||||
|
||||
|
@ -21,11 +21,15 @@
|
||||
// You should have received a copy of the GNU Lesser General Public
|
||||
// License and a copy of the GNU General Public License along with
|
||||
// Eigen. If not, see <http://www.gnu.org/licenses/>.
|
||||
|
||||
#define EIGEN_DONT_VECTORIZE
|
||||
#include "main.h"
|
||||
#include <Eigen/Cholesky>
|
||||
#include <Eigen/LU>
|
||||
|
||||
#ifdef HAS_GSL
|
||||
#include "gsl_helper.h"
|
||||
#endif
|
||||
|
||||
template<typename MatrixType> void cholesky(const MatrixType& m)
|
||||
{
|
||||
/* this test covers the following files:
|
||||
@ -39,38 +43,79 @@ template<typename MatrixType> void cholesky(const MatrixType& m)
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
||||
|
||||
MatrixType a = test_random_matrix<MatrixType>(rows,cols);
|
||||
MatrixType a0 = test_random_matrix<MatrixType>(rows,cols);
|
||||
VectorType vecB = test_random_matrix<VectorType>(rows);
|
||||
MatrixType matB = test_random_matrix<MatrixType>(rows,cols);
|
||||
SquareMatrixType covMat = a * a.adjoint();
|
||||
SquareMatrixType symm = a0 * a0.adjoint();
|
||||
// let's make sure the matrix is not singular or near singular
|
||||
MatrixType a1 = test_random_matrix<MatrixType>(rows,cols);
|
||||
symm += a1 * a1.adjoint();
|
||||
|
||||
#ifdef HAS_GSL
|
||||
if (ei_is_same_type<RealScalar,double>::ret)
|
||||
{
|
||||
typedef GslTraits<Scalar> Gsl;
|
||||
typename Gsl::Matrix gMatA=0, gSymm=0;
|
||||
typename Gsl::Vector gVecB=0, gVecX=0;
|
||||
convert<MatrixType>(symm, gSymm);
|
||||
convert<MatrixType>(symm, gMatA);
|
||||
convert<VectorType>(vecB, gVecB);
|
||||
convert<VectorType>(vecB, gVecX);
|
||||
Gsl::cholesky(gMatA);
|
||||
Gsl::cholesky_solve(gMatA, gVecB, gVecX);
|
||||
VectorType vecX, _vecX, _vecB;
|
||||
convert(gVecX, _vecX);
|
||||
vecX = symm.cholesky().solve(vecB);
|
||||
Gsl::prod(gSymm, gVecX, gVecB);
|
||||
convert(gVecB, _vecB);
|
||||
// test gsl itself !
|
||||
VERIFY_IS_APPROX(vecB, _vecB);
|
||||
VERIFY_IS_APPROX(vecX, _vecX);
|
||||
|
||||
Gsl::free(gMatA);
|
||||
Gsl::free(gSymm);
|
||||
Gsl::free(gVecB);
|
||||
Gsl::free(gVecX);
|
||||
}
|
||||
#endif
|
||||
|
||||
if (rows>1)
|
||||
{
|
||||
CholeskyWithoutSquareRoot<SquareMatrixType> cholnosqrt(covMat);
|
||||
VERIFY_IS_APPROX(covMat, cholnosqrt.matrixL() * cholnosqrt.vectorD().asDiagonal() * cholnosqrt.matrixL().adjoint());
|
||||
// cout << (covMat * cholnosqrt.solve(vecB)).transpose().format(6) << endl;
|
||||
// cout << vecB.transpose().format(6) << endl << "----------" << endl;
|
||||
VERIFY((covMat * cholnosqrt.solve(vecB)).isApprox(vecB, test_precision<RealScalar>()*RealScalar(100))); // FIXME
|
||||
VERIFY((covMat * cholnosqrt.solve(matB)).isApprox(matB, test_precision<RealScalar>()*RealScalar(100))); // FIXME
|
||||
CholeskyWithoutSquareRoot<SquareMatrixType> cholnosqrt(symm);
|
||||
VERIFY(cholnosqrt.isPositiveDefinite());
|
||||
VERIFY_IS_APPROX(symm, cholnosqrt.matrixL() * cholnosqrt.vectorD().asDiagonal() * cholnosqrt.matrixL().adjoint());
|
||||
VERIFY_IS_APPROX(symm * cholnosqrt.solve(vecB), vecB);
|
||||
VERIFY_IS_APPROX(symm * cholnosqrt.solve(matB), matB);
|
||||
}
|
||||
|
||||
Cholesky<SquareMatrixType> chol(covMat);
|
||||
VERIFY_IS_APPROX(covMat, chol.matrixL() * chol.matrixL().adjoint());
|
||||
// cout << (covMat * chol.solve(vecB)).transpose().format(6) << endl;
|
||||
// cout << vecB.transpose().format(6) << endl << "----------" << endl;
|
||||
VERIFY((covMat * chol.solve(vecB)).isApprox(vecB, test_precision<RealScalar>()*RealScalar(100))); // FIXME
|
||||
VERIFY((covMat * chol.solve(matB)).isApprox(matB, test_precision<RealScalar>()*RealScalar(100))); // FIXME
|
||||
{
|
||||
Cholesky<SquareMatrixType> chol(symm);
|
||||
VERIFY(chol.isPositiveDefinite());
|
||||
VERIFY_IS_APPROX(symm, chol.matrixL() * chol.matrixL().adjoint());
|
||||
VERIFY_IS_APPROX(symm * chol.solve(vecB), vecB);
|
||||
VERIFY_IS_APPROX(symm * chol.solve(matB), matB);
|
||||
}
|
||||
|
||||
// test isPositiveDefinite on non definite matrix
|
||||
if (rows>4)
|
||||
{
|
||||
SquareMatrixType symm = a0.block(0,0,rows,cols-4) * a0.block(0,0,rows,cols-4).adjoint();
|
||||
Cholesky<SquareMatrixType> chol(symm);
|
||||
VERIFY(!chol.isPositiveDefinite());
|
||||
CholeskyWithoutSquareRoot<SquareMatrixType> cholnosqrt(symm);
|
||||
VERIFY(!cholnosqrt.isPositiveDefinite());
|
||||
}
|
||||
}
|
||||
|
||||
void test_cholesky()
|
||||
{
|
||||
for(int i = 0; i < g_repeat; i++) {
|
||||
CALL_SUBTEST( cholesky(Matrix<float,1,1>()) );
|
||||
CALL_SUBTEST( cholesky(Matrix<float,2,2>()) );
|
||||
// CALL_SUBTEST( cholesky(Matrix3f()) );
|
||||
// CALL_SUBTEST( cholesky(Matrix4d()) );
|
||||
// CALL_SUBTEST( cholesky(MatrixXcd(7,7)) );
|
||||
// CALL_SUBTEST( cholesky(MatrixXf(19,19)) );
|
||||
// CALL_SUBTEST( cholesky(MatrixXd(33,33)) );
|
||||
CALL_SUBTEST( cholesky(Matrix<double,1,1>()) );
|
||||
CALL_SUBTEST( cholesky(Matrix2d()) );
|
||||
CALL_SUBTEST( cholesky(Matrix3f()) );
|
||||
CALL_SUBTEST( cholesky(Matrix4d()) );
|
||||
CALL_SUBTEST( cholesky(MatrixXcd(7,7)) );
|
||||
CALL_SUBTEST( cholesky(MatrixXf(17,17)) );
|
||||
CALL_SUBTEST( cholesky(MatrixXd(33,33)) );
|
||||
}
|
||||
}
|
||||
|
@ -25,6 +25,10 @@
|
||||
#include "main.h"
|
||||
#include <Eigen/QR>
|
||||
|
||||
#ifdef HAS_GSL
|
||||
#include "gsl_helper.h"
|
||||
#endif
|
||||
|
||||
template<typename MatrixType> void eigensolver(const MatrixType& m)
|
||||
{
|
||||
/* this test covers the following files:
|
||||
@ -33,19 +37,76 @@ template<typename MatrixType> void eigensolver(const MatrixType& m)
|
||||
int rows = m.rows();
|
||||
int cols = m.cols();
|
||||
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
||||
typedef Matrix<RealScalar, MatrixType::RowsAtCompileTime, 1> RealVectorType;
|
||||
typedef typename std::complex<typename NumTraits<typename MatrixType::Scalar>::Real> Complex;
|
||||
|
||||
MatrixType a = MatrixType::Random(rows,cols);
|
||||
MatrixType symmA = a.adjoint() * a;
|
||||
RealScalar largerEps = 10*test_precision<RealScalar>();
|
||||
|
||||
MatrixType a = test_random_matrix<MatrixType>(rows,cols);
|
||||
MatrixType a1 = test_random_matrix<MatrixType>(rows,cols);
|
||||
MatrixType symmA = a.adjoint() * a + a1.adjoint() * a1;
|
||||
|
||||
MatrixType b = test_random_matrix<MatrixType>(rows,cols);
|
||||
MatrixType b1 = test_random_matrix<MatrixType>(rows,cols);
|
||||
MatrixType symmB = b.adjoint() * b + b1.adjoint() * b1;
|
||||
|
||||
SelfAdjointEigenSolver<MatrixType> eiSymm(symmA);
|
||||
VERIFY_IS_APPROX(symmA * eiSymm.eigenvectors(), (eiSymm.eigenvectors() * eiSymm.eigenvalues().asDiagonal().eval()));
|
||||
// generalized eigen pb
|
||||
SelfAdjointEigenSolver<MatrixType> eiSymmGen(symmA, symmB);
|
||||
|
||||
#ifdef HAS_GSL
|
||||
if (ei_is_same_type<RealScalar,double>::ret)
|
||||
{
|
||||
typedef GslTraits<Scalar> Gsl;
|
||||
typename Gsl::Matrix gEvec=0, gSymmA=0, gSymmB=0;
|
||||
typename GslTraits<RealScalar>::Vector gEval=0;
|
||||
RealVectorType _eval;
|
||||
MatrixType _evec;
|
||||
convert<MatrixType>(symmA, gSymmA);
|
||||
convert<MatrixType>(symmB, gSymmB);
|
||||
convert<MatrixType>(symmA, gEvec);
|
||||
gEval = GslTraits<RealScalar>::createVector(rows);
|
||||
|
||||
Gsl::eigen_symm(gSymmA, gEval, gEvec);
|
||||
convert(gEval, _eval);
|
||||
convert(gEvec, _evec);
|
||||
|
||||
// test gsl itself !
|
||||
VERIFY((symmA * _evec).isApprox(_evec * _eval.asDiagonal().eval(), largerEps));
|
||||
|
||||
// compare with eigen
|
||||
VERIFY_IS_APPROX(_eval, eiSymm.eigenvalues());
|
||||
VERIFY_IS_APPROX(_evec.cwise().abs(), eiSymm.eigenvectors().cwise().abs());
|
||||
|
||||
// generalized pb
|
||||
Gsl::eigen_symm_gen(gSymmA, gSymmB, gEval, gEvec);
|
||||
convert(gEval, _eval);
|
||||
convert(gEvec, _evec);
|
||||
// test GSL itself:
|
||||
VERIFY((symmA * _evec).isApprox(symmB * (_evec * _eval.asDiagonal().eval()), largerEps));
|
||||
|
||||
// compare with eigen
|
||||
// std::cerr << _eval.transpose() << "\n" << eiSymmGen.eigenvalues().transpose() << "\n\n";
|
||||
// std::cerr << _evec.format(6) << "\n\n" << eiSymmGen.eigenvectors().format(6) << "\n\n\n";
|
||||
VERIFY_IS_APPROX(_eval, eiSymmGen.eigenvalues());
|
||||
VERIFY_IS_APPROX(_evec.cwise().abs(), eiSymmGen.eigenvectors().cwise().abs());
|
||||
|
||||
Gsl::free(gSymmA);
|
||||
Gsl::free(gSymmB);
|
||||
GslTraits<RealScalar>::free(gEval);
|
||||
Gsl::free(gEvec);
|
||||
}
|
||||
#endif
|
||||
|
||||
VERIFY((symmA * eiSymm.eigenvectors()).isApprox(
|
||||
eiSymm.eigenvectors() * eiSymm.eigenvalues().asDiagonal().eval(), largerEps));
|
||||
|
||||
// generalized eigen problem Ax = lBx
|
||||
MatrixType b = MatrixType::Random(rows,cols);
|
||||
MatrixType symmB = b.adjoint() * b;
|
||||
eiSymm.compute(symmA,symmB);
|
||||
VERIFY_IS_APPROX(symmA * eiSymm.eigenvectors(), symmB * (eiSymm.eigenvectors() * eiSymm.eigenvalues().asDiagonal().eval()));
|
||||
VERIFY((symmA * eiSymmGen.eigenvectors()).isApprox(
|
||||
symmB * (eiSymmGen.eigenvectors() * eiSymmGen.eigenvalues().asDiagonal().eval()), largerEps));
|
||||
|
||||
// EigenSolver<MatrixType> eiNotSymmButSymm(covMat);
|
||||
// VERIFY_IS_APPROX((covMat.template cast<Complex>()) * (eiNotSymmButSymm.eigenvectors().template cast<Complex>()),
|
||||
@ -59,12 +120,12 @@ template<typename MatrixType> void eigensolver(const MatrixType& m)
|
||||
|
||||
void test_eigensolver()
|
||||
{
|
||||
for(int i = 0; i < 1; i++) {
|
||||
for(int i = 0; i < g_repeat; i++) {
|
||||
// very important to test a 3x3 matrix since we provide a special path for it
|
||||
CALL_SUBTEST( eigensolver(Matrix3f()) );
|
||||
CALL_SUBTEST( eigensolver(Matrix4d()) );
|
||||
CALL_SUBTEST( eigensolver(MatrixXf(7,7)) );
|
||||
CALL_SUBTEST( eigensolver(MatrixXcd(6,6)) );
|
||||
CALL_SUBTEST( eigensolver(MatrixXcf(3,3)) );
|
||||
CALL_SUBTEST( eigensolver(MatrixXcd(5,5)) );
|
||||
CALL_SUBTEST( eigensolver(MatrixXd(19,19)) );
|
||||
}
|
||||
}
|
||||
|
@ -69,8 +69,8 @@ template<typename Scalar> void geometry(void)
|
||||
VERIFY_IS_APPROX(q1 * v2, q1.toRotationMatrix() * v2);
|
||||
VERIFY_IS_APPROX(q1 * q2 * v2,
|
||||
q1.toRotationMatrix() * q2.toRotationMatrix() * v2);
|
||||
VERIFY_IS_NOT_APPROX(q2 * q1 * v2,
|
||||
q1.toRotationMatrix() * q2.toRotationMatrix() * v2);
|
||||
VERIFY( !(q2 * q1 * v2).isApprox(
|
||||
q1.toRotationMatrix() * q2.toRotationMatrix() * v2));
|
||||
q2 = q1.toRotationMatrix();
|
||||
VERIFY_IS_APPROX(q1*v1,q2*v1);
|
||||
|
||||
@ -126,7 +126,7 @@ template<typename Scalar> void geometry(void)
|
||||
t1.prescale(v0);
|
||||
|
||||
VERIFY_IS_APPROX( (t0 * Vector3(1,0,0)).norm(), v0.x());
|
||||
VERIFY_IS_NOT_APPROX((t1 * Vector3(1,0,0)).norm(), v0.x());
|
||||
VERIFY(!ei_isApprox((t1 * Vector3(1,0,0)).norm(), v0.x()));
|
||||
|
||||
t0.setIdentity();
|
||||
t1.setIdentity();
|
||||
@ -138,7 +138,7 @@ template<typename Scalar> void geometry(void)
|
||||
t1.prescale(v1.cwise().inverse());
|
||||
t1.translate(-v0);
|
||||
|
||||
VERIFY((t0.matrix() * t1.matrix()).isIdentity());
|
||||
VERIFY((t0.matrix() * t1.matrix()).isIdentity(test_precision<Scalar>()));
|
||||
|
||||
t1.fromPositionOrientationScale(v0, q1, v1);
|
||||
VERIFY_IS_APPROX(t1.matrix(), t0.matrix());
|
||||
@ -147,6 +147,8 @@ template<typename Scalar> void geometry(void)
|
||||
Transform2 t20, t21;
|
||||
Vector2 v20 = test_random_matrix<Vector2>();
|
||||
Vector2 v21 = test_random_matrix<Vector2>();
|
||||
for (int k=0; k<2; ++k)
|
||||
if (ei_abs(v21[k])<1e-3) v21[k] = 1e-3;
|
||||
t21.setIdentity();
|
||||
t21.linear() = Rotation2D<Scalar>(a).toRotationMatrix();
|
||||
VERIFY_IS_APPROX(t20.fromPositionOrientationScale(v20,a,v21).matrix(),
|
||||
@ -154,7 +156,8 @@ template<typename Scalar> void geometry(void)
|
||||
|
||||
t21.setIdentity();
|
||||
t21.linear() = Rotation2D<Scalar>(-a).toRotationMatrix();
|
||||
VERIFY( (t20.fromPositionOrientationScale(v20,a,v21) * (t21.prescale(v21.cwise().inverse()).translate(-v20))).isIdentity() );
|
||||
VERIFY( (t20.fromPositionOrientationScale(v20,a,v21)
|
||||
* (t21.prescale(v21.cwise().inverse()).translate(-v20))).isIdentity(test_precision<Scalar>()) );
|
||||
}
|
||||
|
||||
void test_geometry()
|
||||
|
190
test/gsl_helper.h
Normal file
190
test/gsl_helper.h
Normal file
@ -0,0 +1,190 @@
|
||||
// This file is part of Eigen, a lightweight C++ template library
|
||||
// for linear algebra. Eigen itself is part of the KDE project.
|
||||
//
|
||||
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
|
||||
//
|
||||
// Eigen is free software; you can redistribute it and/or
|
||||
// modify it under the terms of the GNU Lesser General Public
|
||||
// License as published by the Free Software Foundation; either
|
||||
// version 3 of the License, or (at your option) any later version.
|
||||
//
|
||||
// Alternatively, you can redistribute it and/or
|
||||
// modify it under the terms of the GNU General Public License as
|
||||
// published by the Free Software Foundation; either version 2 of
|
||||
// the License, or (at your option) any later version.
|
||||
//
|
||||
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
|
||||
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
|
||||
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
|
||||
// GNU General Public License for more details.
|
||||
//
|
||||
// You should have received a copy of the GNU Lesser General Public
|
||||
// License and a copy of the GNU General Public License along with
|
||||
// Eigen. If not, see <http://www.gnu.org/licenses/>.
|
||||
|
||||
#ifndef EIGEN_GSL_HELPER
|
||||
#define EIGEN_GSL_HELPER
|
||||
|
||||
#include <Eigen/Core>
|
||||
|
||||
#include <gsl/gsl_blas.h>
|
||||
#include <gsl/gsl_multifit.h>
|
||||
#include <gsl/gsl_eigen.h>
|
||||
#include <gsl/gsl_linalg.h>
|
||||
#include <gsl/gsl_complex.h>
|
||||
#include <gsl/gsl_complex_math.h>
|
||||
|
||||
namespace Eigen {
|
||||
|
||||
template<typename Scalar, bool IsComplex = NumTraits<Scalar>::IsComplex> struct GslTraits
|
||||
{
|
||||
typedef gsl_matrix* Matrix;
|
||||
typedef gsl_vector* Vector;
|
||||
static Matrix createMatrix(int rows, int cols) { return gsl_matrix_alloc(rows,cols); }
|
||||
static Vector createVector(int size) { return gsl_vector_alloc(size); }
|
||||
static void free(Matrix& m) { gsl_matrix_free(m); m=0; }
|
||||
static void free(Vector& m) { gsl_vector_free(m); m=0; }
|
||||
static void prod(const Matrix& m, const Vector& v, Vector& x) { gsl_blas_dgemv(CblasNoTrans,1,m,v,0,x); }
|
||||
static void cholesky(Matrix& m) { gsl_linalg_cholesky_decomp(m); }
|
||||
static void cholesky_solve(const Matrix& m, const Vector& b, Vector& x) { gsl_linalg_cholesky_solve(m,b,x); }
|
||||
static void eigen_symm(const Matrix& m, Vector& eval, Matrix& evec)
|
||||
{
|
||||
gsl_eigen_symmv_workspace * w = gsl_eigen_symmv_alloc(m->size1);
|
||||
Matrix a = createMatrix(m->size1, m->size2);
|
||||
gsl_matrix_memcpy(a, m);
|
||||
gsl_eigen_symmv(a,eval,evec,w);
|
||||
gsl_eigen_symmv_sort(eval, evec, GSL_EIGEN_SORT_VAL_ASC);
|
||||
gsl_eigen_symmv_free(w);
|
||||
free(a);
|
||||
}
|
||||
static void eigen_symm_gen(const Matrix& m, const Matrix& _b, Vector& eval, Matrix& evec)
|
||||
{
|
||||
gsl_eigen_gensymmv_workspace * w = gsl_eigen_gensymmv_alloc(m->size1);
|
||||
Matrix a = createMatrix(m->size1, m->size2);
|
||||
Matrix b = createMatrix(_b->size1, _b->size2);
|
||||
gsl_matrix_memcpy(a, m);
|
||||
gsl_matrix_memcpy(b, _b);
|
||||
gsl_eigen_gensymmv(a,b,eval,evec,w);
|
||||
gsl_eigen_symmv_sort(eval, evec, GSL_EIGEN_SORT_VAL_ASC);
|
||||
gsl_eigen_gensymmv_free(w);
|
||||
free(a);
|
||||
}
|
||||
};
|
||||
|
||||
template<typename Scalar> struct GslTraits<Scalar,true>
|
||||
{
|
||||
typedef gsl_matrix_complex* Matrix;
|
||||
typedef gsl_vector_complex* Vector;
|
||||
static Matrix createMatrix(int rows, int cols) { return gsl_matrix_complex_alloc(rows,cols); }
|
||||
static Vector createVector(int size) { return gsl_vector_complex_alloc(size); }
|
||||
static void free(Matrix& m) { gsl_matrix_complex_free(m); m=0; }
|
||||
static void free(Vector& m) { gsl_vector_complex_free(m); m=0; }
|
||||
static void cholesky(Matrix& m) { gsl_linalg_complex_cholesky_decomp(m); }
|
||||
static void cholesky_solve(const Matrix& m, const Vector& b, Vector& x) { gsl_linalg_complex_cholesky_solve(m,b,x); }
|
||||
static void prod(const Matrix& m, const Vector& v, Vector& x)
|
||||
{ gsl_blas_zgemv(CblasNoTrans,gsl_complex_rect(1,0),m,v,gsl_complex_rect(0,0),x); }
|
||||
static void eigen_symm(const Matrix& m, gsl_vector* &eval, Matrix& evec)
|
||||
{
|
||||
gsl_eigen_hermv_workspace * w = gsl_eigen_hermv_alloc(m->size1);
|
||||
Matrix a = createMatrix(m->size1, m->size2);
|
||||
gsl_matrix_complex_memcpy(a, m);
|
||||
gsl_eigen_hermv(a,eval,evec,w);
|
||||
gsl_eigen_hermv_sort(eval, evec, GSL_EIGEN_SORT_VAL_ASC);
|
||||
gsl_eigen_hermv_free(w);
|
||||
free(a);
|
||||
}
|
||||
static void eigen_symm_gen(const Matrix& m, const Matrix& _b, gsl_vector* &eval, Matrix& evec)
|
||||
{
|
||||
gsl_eigen_genhermv_workspace * w = gsl_eigen_genhermv_alloc(m->size1);
|
||||
Matrix a = createMatrix(m->size1, m->size2);
|
||||
Matrix b = createMatrix(_b->size1, _b->size2);
|
||||
gsl_matrix_complex_memcpy(a, m);
|
||||
gsl_matrix_complex_memcpy(b, _b);
|
||||
gsl_eigen_genhermv(a,b,eval,evec,w);
|
||||
gsl_eigen_hermv_sort(eval, evec, GSL_EIGEN_SORT_VAL_ASC);
|
||||
gsl_eigen_genhermv_free(w);
|
||||
free(a);
|
||||
}
|
||||
};
|
||||
|
||||
template<typename MatrixType>
|
||||
void convert(const MatrixType& m, gsl_matrix* &res)
|
||||
{
|
||||
// if (res)
|
||||
// gsl_matrix_free(res);
|
||||
res = gsl_matrix_alloc(m.rows(), m.cols());
|
||||
for (int i=0 ; i<m.rows() ; ++i)
|
||||
for (int j=0 ; j<m.cols(); ++j)
|
||||
gsl_matrix_set(res, i, j, m(i,j));
|
||||
}
|
||||
|
||||
template<typename MatrixType>
|
||||
void convert(const gsl_matrix* m, MatrixType& res)
|
||||
{
|
||||
res.resize(int(m->size1), int(m->size2));
|
||||
for (int i=0 ; i<res.rows() ; ++i)
|
||||
for (int j=0 ; j<res.cols(); ++j)
|
||||
res(i,j) = gsl_matrix_get(m,i,j);
|
||||
}
|
||||
|
||||
template<typename VectorType>
|
||||
void convert(const VectorType& m, gsl_vector* &res)
|
||||
{
|
||||
if (res) gsl_vector_free(res);
|
||||
res = gsl_vector_alloc(m.size());
|
||||
for (int i=0 ; i<m.size() ; ++i)
|
||||
gsl_vector_set(res, i, m[i]);
|
||||
}
|
||||
|
||||
template<typename VectorType>
|
||||
void convert(const gsl_vector* m, VectorType& res)
|
||||
{
|
||||
res.resize (m->size);
|
||||
for (int i=0 ; i<res.rows() ; ++i)
|
||||
res[i] = gsl_vector_get(m, i);
|
||||
}
|
||||
|
||||
template<typename MatrixType>
|
||||
void convert(const MatrixType& m, gsl_matrix_complex* &res)
|
||||
{
|
||||
res = gsl_matrix_complex_alloc(m.rows(), m.cols());
|
||||
for (int i=0 ; i<m.rows() ; ++i)
|
||||
for (int j=0 ; j<m.cols(); ++j)
|
||||
{
|
||||
gsl_matrix_complex_set(res, i, j,
|
||||
gsl_complex_rect(m(i,j).real(), m(i,j).imag()));
|
||||
}
|
||||
}
|
||||
|
||||
template<typename MatrixType>
|
||||
void convert(const gsl_matrix_complex* m, MatrixType& res)
|
||||
{
|
||||
res.resize(int(m->size1), int(m->size2));
|
||||
for (int i=0 ; i<res.rows() ; ++i)
|
||||
for (int j=0 ; j<res.cols(); ++j)
|
||||
res(i,j) = typename MatrixType::Scalar(
|
||||
GSL_REAL(gsl_matrix_complex_get(m,i,j)),
|
||||
GSL_IMAG(gsl_matrix_complex_get(m,i,j)));
|
||||
}
|
||||
|
||||
template<typename VectorType>
|
||||
void convert(const VectorType& m, gsl_vector_complex* &res)
|
||||
{
|
||||
res = gsl_vector_complex_alloc(m.size());
|
||||
for (int i=0 ; i<m.size() ; ++i)
|
||||
gsl_vector_complex_set(res, i, gsl_complex_rect(m[i].real(), m[i].imag()));
|
||||
}
|
||||
|
||||
template<typename VectorType>
|
||||
void convert(const gsl_vector_complex* m, VectorType& res)
|
||||
{
|
||||
res.resize(m->size);
|
||||
for (int i=0 ; i<res.rows() ; ++i)
|
||||
res[i] = typename VectorType::Scalar(
|
||||
GSL_REAL(gsl_vector_complex_get(m, i)),
|
||||
GSL_IMAG(gsl_vector_complex_get(m, i)));
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
#endif // EIGEN_GSL_HELPER
|
@ -35,13 +35,21 @@ template<typename MatrixType> void inverse(const MatrixType& m)
|
||||
int cols = m.cols();
|
||||
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> VectorType;
|
||||
|
||||
MatrixType m1 = test_random_matrix<MatrixType>(rows, cols),
|
||||
m2 = test_random_matrix<MatrixType>(rows, cols),
|
||||
m2(rows, cols),
|
||||
mzero = MatrixType::Zero(rows, cols),
|
||||
identity = MatrixType::Identity(rows, rows);
|
||||
|
||||
if (ei_is_same_type<RealScalar,float>::ret)
|
||||
{
|
||||
// let's build a more stable to inverse matrix
|
||||
MatrixType a = test_random_matrix<MatrixType>(rows,cols);
|
||||
m1 += m1 * m1.adjoint() + a * a.adjoint();
|
||||
}
|
||||
|
||||
m2 = m1.inverse();
|
||||
VERIFY_IS_APPROX(m1, m2.inverse() );
|
||||
|
||||
|
@ -41,15 +41,10 @@ template<typename MatrixType> void linearStructure(const MatrixType& m)
|
||||
MatrixType m1 = test_random_matrix<MatrixType>(rows, cols),
|
||||
m2 = test_random_matrix<MatrixType>(rows, cols),
|
||||
m3(rows, cols),
|
||||
mzero = MatrixType::Zero(rows, cols),
|
||||
identity = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
::Identity(rows, rows),
|
||||
square = test_random_matrix<Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> >(rows, rows);
|
||||
VectorType v1 = test_random_matrix<VectorType>(rows),
|
||||
v2 = test_random_matrix<VectorType>(rows),
|
||||
vzero = VectorType::Zero(rows);
|
||||
mzero = MatrixType::Zero(rows, cols);
|
||||
|
||||
Scalar s1 = test_random<Scalar>();
|
||||
while (ei_abs(s1)<1e-3) s1 = test_random<Scalar>();
|
||||
|
||||
int r = ei_random<int>(0, rows-1),
|
||||
c = ei_random<int>(0, cols-1);
|
||||
@ -94,6 +89,7 @@ void test_linearstructure()
|
||||
for(int i = 0; i < g_repeat; i++) {
|
||||
CALL_SUBTEST( linearStructure(Matrix<float, 1, 1>()) );
|
||||
CALL_SUBTEST( linearStructure(Matrix2f()) );
|
||||
CALL_SUBTEST( linearStructure(Vector3d()) );
|
||||
CALL_SUBTEST( linearStructure(Matrix4d()) );
|
||||
CALL_SUBTEST( linearStructure(MatrixXcf(3, 3)) );
|
||||
CALL_SUBTEST( linearStructure(MatrixXf(8, 12)) );
|
||||
|
12
test/lu.cpp
12
test/lu.cpp
@ -51,7 +51,8 @@ template<typename MatrixType> void lu_non_invertible()
|
||||
/* this test covers the following files:
|
||||
LU.h
|
||||
*/
|
||||
int rows = ei_random<int>(10,200), cols = ei_random<int>(10,200), cols2 = ei_random<int>(10,200);
|
||||
// NOTE lu.dimensionOfKernel() fails most of the time for rows or cols smaller that 11
|
||||
int rows = ei_random<int>(11,200), cols = ei_random<int>(11,200), cols2 = ei_random<int>(11,200);
|
||||
int rank = ei_random<int>(1, std::min(rows, cols)-1);
|
||||
|
||||
MatrixType m1(rows, cols), m2(cols, cols2), m3(rows, cols2), k(1,1);
|
||||
@ -91,6 +92,13 @@ template<typename MatrixType> void lu_invertible()
|
||||
MatrixType m1(size, size), m2(size, size), m3(size, size);
|
||||
m1 = test_random_matrix<MatrixType>(size,size);
|
||||
|
||||
if (ei_is_same_type<RealScalar,float>::ret)
|
||||
{
|
||||
// let's build a matrix more stable to inverse
|
||||
MatrixType a = test_random_matrix<MatrixType>(size,size*2);
|
||||
m1 += a * a.adjoint();
|
||||
}
|
||||
|
||||
LU<MatrixType> lu(m1);
|
||||
VERIFY(0 == lu.dimensionOfKernel());
|
||||
VERIFY(size == lu.rank());
|
||||
@ -99,7 +107,7 @@ template<typename MatrixType> void lu_invertible()
|
||||
VERIFY(lu.isInvertible());
|
||||
m3 = test_random_matrix<MatrixType>(size,size);
|
||||
lu.solve(m3, &m2);
|
||||
VERIFY(m3.isApprox(m1*m2, test_precision<RealScalar>()*RealScalar(100))); // FIXME
|
||||
VERIFY_IS_APPROX(m3, m1*m2);
|
||||
VERIFY_IS_APPROX(m2, lu.inverse()*m3);
|
||||
m3 = test_random_matrix<MatrixType>(size,size);
|
||||
VERIFY(lu.solve(m3, &m2));
|
||||
|
@ -29,6 +29,7 @@ void test_product_small()
|
||||
for(int i = 0; i < g_repeat; i++) {
|
||||
CALL_SUBTEST( product(Matrix<float, 3, 2>()) );
|
||||
CALL_SUBTEST( product(Matrix<int, 3, 5>()) );
|
||||
CALL_SUBTEST( product(Matrix3d()) );
|
||||
CALL_SUBTEST( product(Matrix4d()) );
|
||||
CALL_SUBTEST( product(Matrix4f()) );
|
||||
}
|
||||
|
@ -10,7 +10,7 @@ cyan='\E[36m'
|
||||
white='\E[37m'
|
||||
|
||||
if make test_$1 > /dev/null 2> .runtest.log ; then
|
||||
if ! ./test_$1 > /dev/null 2> .runtest.log ; then
|
||||
if ! ./test_$1 r20 > /dev/null 2> .runtest.log ; then
|
||||
echo -e $red Test $1 failed: $black
|
||||
echo -e $blue
|
||||
cat .runtest.log
|
||||
|
19
test/svd.cpp
19
test/svd.cpp
@ -34,11 +34,16 @@ template<typename MatrixType> void svd(const MatrixType& m)
|
||||
int cols = m.cols();
|
||||
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
MatrixType a = MatrixType::Random(rows,cols);
|
||||
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||
MatrixType a = test_random_matrix<MatrixType>(rows,cols);
|
||||
Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> b =
|
||||
Matrix<Scalar, MatrixType::RowsAtCompileTime, 1>::Random(rows,1);
|
||||
test_random_matrix<Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> >(rows,1);
|
||||
Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> x(cols,1), x2(cols,1);
|
||||
|
||||
RealScalar largerEps = test_precision<RealScalar>();
|
||||
if (ei_is_same_type<RealScalar,float>::ret)
|
||||
largerEps = 1e-3f;
|
||||
|
||||
SVD<MatrixType> svd(a);
|
||||
MatrixType sigma = MatrixType::Zero(rows,cols);
|
||||
MatrixType matU = MatrixType::Zero(rows,rows);
|
||||
@ -49,8 +54,14 @@ template<typename MatrixType> void svd(const MatrixType& m)
|
||||
|
||||
if (rows==cols)
|
||||
{
|
||||
if (ei_is_same_type<RealScalar,float>::ret)
|
||||
{
|
||||
MatrixType a1 = test_random_matrix<MatrixType>(rows,cols);
|
||||
a += a * a.adjoint() + a1 * a1.adjoint();
|
||||
}
|
||||
SVD<MatrixType> svd(a);
|
||||
svd.solve(b, &x);
|
||||
VERIFY_IS_APPROX(a * x, b);
|
||||
VERIFY_IS_APPROX(a * x,b);
|
||||
}
|
||||
}
|
||||
|
||||
@ -60,7 +71,7 @@ void test_svd()
|
||||
CALL_SUBTEST( svd(Matrix3f()) );
|
||||
CALL_SUBTEST( svd(Matrix4d()) );
|
||||
CALL_SUBTEST( svd(MatrixXf(7,7)) );
|
||||
CALL_SUBTEST( svd(MatrixXf(14,7)) );
|
||||
CALL_SUBTEST( svd(MatrixXd(14,7)) );
|
||||
// complex are not implemented yet
|
||||
// CALL_SUBTEST( svd(MatrixXcd(6,6)) );
|
||||
// CALL_SUBTEST( svd(MatrixXcf(3,3)) );
|
||||
|
@ -30,12 +30,15 @@ template<typename MatrixType> void triangular(const MatrixType& m)
|
||||
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
||||
|
||||
RealScalar largerEps = 10*test_precision<RealScalar>();
|
||||
|
||||
int rows = m.rows();
|
||||
int cols = m.cols();
|
||||
|
||||
MatrixType m1 = MatrixType::Random(rows, cols),
|
||||
m2 = MatrixType::Random(rows, cols),
|
||||
MatrixType m1 = test_random_matrix<MatrixType>(rows, cols),
|
||||
m2 = test_random_matrix<MatrixType>(rows, cols),
|
||||
m3(rows, cols),
|
||||
m4(rows, cols),
|
||||
r1(rows, cols),
|
||||
r2(rows, cols),
|
||||
mzero = MatrixType::Zero(rows, cols),
|
||||
@ -44,8 +47,8 @@ template<typename MatrixType> void triangular(const MatrixType& m)
|
||||
::Identity(rows, rows),
|
||||
square = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
::Random(rows, rows);
|
||||
VectorType v1 = VectorType::Random(rows),
|
||||
v2 = VectorType::Random(rows),
|
||||
VectorType v1 = test_random_matrix<VectorType>(rows),
|
||||
v2 = test_random_matrix<VectorType>(rows),
|
||||
vzero = VectorType::Zero(rows);
|
||||
|
||||
MatrixType m1up = m1.template part<Eigen::Upper>();
|
||||
@ -78,17 +81,34 @@ template<typename MatrixType> void triangular(const MatrixType& m)
|
||||
m1.template part<Eigen::Lower>() = (m2.transpose() * m2).lazy();
|
||||
VERIFY_IS_APPROX(m3.template part<Eigen::Lower>(), m1);
|
||||
|
||||
m1 = test_random_matrix<MatrixType>(rows, cols);
|
||||
for (int i=0; i<rows; ++i)
|
||||
while (ei_abs2(m1(i,i))<1e-3) m1(i,i) = test_random<Scalar>();
|
||||
|
||||
Transpose<MatrixType> trm4(m4);
|
||||
// test back and forward subsitution
|
||||
m3 = m1.template part<Eigen::Lower>();
|
||||
VERIFY(m3.template marked<Eigen::Lower>().solveTriangular(m3).cwise().abs().isIdentity(test_precision<RealScalar>()));
|
||||
VERIFY(m3.transpose().template marked<Eigen::Upper>()
|
||||
.solveTriangular(m3.transpose()).cwise().abs().isIdentity(test_precision<RealScalar>()));
|
||||
// check M * inv(L) using in place API
|
||||
m4 = m3;
|
||||
m3.transpose().template marked<Eigen::Upper>().solveTriangularInPlace(trm4);
|
||||
VERIFY(m4.cwise().abs().isIdentity(test_precision<RealScalar>()));
|
||||
|
||||
m3 = m1.template part<Eigen::Upper>();
|
||||
VERIFY(m3.template marked<Eigen::Upper>().solveTriangular(m3).cwise().abs().isIdentity(test_precision<RealScalar>()));
|
||||
VERIFY(m3.transpose().template marked<Eigen::Lower>()
|
||||
.solveTriangular(m3.transpose()).cwise().abs().isIdentity(test_precision<RealScalar>()));
|
||||
// check M * inv(U) using in place API
|
||||
m4 = m3;
|
||||
m3.transpose().template marked<Eigen::Lower>().solveTriangularInPlace(trm4);
|
||||
VERIFY(m4.cwise().abs().isIdentity(test_precision<RealScalar>()));
|
||||
|
||||
// FIXME these tests failed due to numerical issues
|
||||
// m1 = MatrixType::Random(rows, cols);
|
||||
// VERIFY_IS_APPROX(m1.template part<Eigen::Upper>().eval() * (m1.template part<Eigen::Upper>().solveTriangular(m2)), m2);
|
||||
// VERIFY_IS_APPROX(m1.template part<Eigen::Lower>().eval() * (m1.template part<Eigen::Lower>().solveTriangular(m2)), m2);
|
||||
m3 = m1.template part<Eigen::Upper>();
|
||||
VERIFY(m2.isApprox(m3 * (m3.template marked<Eigen::Upper>().solveTriangular(m2)), largerEps));
|
||||
m3 = m1.template part<Eigen::Lower>();
|
||||
VERIFY(m2.isApprox(m3 * (m3.template marked<Eigen::Lower>().solveTriangular(m2)), largerEps));
|
||||
|
||||
VERIFY((m1.template part<Eigen::Upper>() * m2.template part<Eigen::Upper>()).isUpper());
|
||||
|
||||
@ -102,6 +122,6 @@ void test_triangular()
|
||||
CALL_SUBTEST( triangular(Matrix3d()) );
|
||||
CALL_SUBTEST( triangular(MatrixXcf(4, 4)) );
|
||||
CALL_SUBTEST( triangular(Matrix<std::complex<float>,8, 8>()) );
|
||||
CALL_SUBTEST( triangular(MatrixXf(85,85)) );
|
||||
CALL_SUBTEST( triangular(MatrixXd(17,17)) );
|
||||
}
|
||||
}
|
||||
|
Loading…
x
Reference in New Issue
Block a user