// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2021 Kolja Brix // // 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/. #include "main.h" #include template void check_generateRandomUnitaryMatrix(const Index dim) { const MatrixType Q = generateRandomUnitaryMatrix(dim); // validate dimensions VERIFY_IS_EQUAL(Q.rows(), dim); VERIFY_IS_EQUAL(Q.cols(), dim); VERIFY_IS_UNITARY(Q); } template void check_setupRandomSvs(const Index dim, const RealScalarType max) { const VectorType v = setupRandomSvs(dim, max); // validate dimensions VERIFY_IS_EQUAL(v.size(), dim); // check entries for(Index i = 0; i < v.size(); ++i) VERIFY_GE(v(i), 0); for(Index i = 0; i < v.size()-1; ++i) VERIFY_GE(v(i), v(i+1)); } template void check_setupRangeSvs(const Index dim, const RealScalarType min, const RealScalarType max) { const VectorType v = setupRangeSvs(dim, min, max); // validate dimensions VERIFY_IS_EQUAL(v.size(), dim); // check entries if(dim == 1) { VERIFY_IS_APPROX(v(0), min); } else { VERIFY_IS_APPROX(v(0), max); VERIFY_IS_APPROX(v(dim-1), min); } for(Index i = 0; i < v.size()-1; ++i) VERIFY_GE(v(i), v(i+1)); } template void check_generateRandomMatrixSvs(const Index rows, const Index cols, const Index diag_size, const RealScalar min_svs, const RealScalar max_svs) { RealVectorType svs = setupRangeSvs(diag_size, min_svs, max_svs); MatrixType M; generateRandomMatrixSvs(svs, rows, cols, M); // validate dimensions VERIFY_IS_EQUAL(M.rows(), rows); VERIFY_IS_EQUAL(M.cols(), cols); VERIFY_IS_EQUAL(svs.size(), diag_size); // validate singular values Eigen::JacobiSVD SVD(M); VERIFY_IS_APPROX(svs, SVD.singularValues()); } template void check_random_matrix(const MatrixType &m) { enum { Rows = MatrixType::RowsAtCompileTime, Cols = MatrixType::ColsAtCompileTime, DiagSize = EIGEN_SIZE_MIN_PREFER_DYNAMIC(Rows, Cols) }; typedef typename MatrixType::Scalar Scalar; typedef typename NumTraits::Real RealScalar; typedef Matrix RealVectorType; const Index rows = m.rows(), cols = m.cols(); const Index diag_size = (std::min)(rows, cols); const RealScalar min_svs = 1.0, max_svs = 1000.0; // check generation of unitary random matrices typedef Matrix MatrixAType; typedef Matrix MatrixBType; check_generateRandomUnitaryMatrix(rows); check_generateRandomUnitaryMatrix(cols); // test generators for singular values check_setupRandomSvs(diag_size, max_svs); check_setupRangeSvs(diag_size, min_svs, max_svs); // check generation of random matrices check_generateRandomMatrixSvs(rows, cols, diag_size, min_svs, max_svs); } EIGEN_DECLARE_TEST(random_matrix) { for(int i = 0; i < g_repeat; i++) { CALL_SUBTEST_1(check_random_matrix(Matrix())); CALL_SUBTEST_2(check_random_matrix(Matrix())); CALL_SUBTEST_3(check_random_matrix(Matrix())); CALL_SUBTEST_4(check_random_matrix(Matrix())); CALL_SUBTEST_5(check_random_matrix(Matrix())); CALL_SUBTEST_6(check_random_matrix(Matrix())); CALL_SUBTEST_7(check_random_matrix(Matrix())); CALL_SUBTEST_8(check_random_matrix(Matrix())); CALL_SUBTEST_9(check_random_matrix(Matrix, 12, 12>())); CALL_SUBTEST_10(check_random_matrix(Matrix, 7, 14>())); CALL_SUBTEST_11(check_random_matrix(Matrix, 15, 11>())); CALL_SUBTEST_12(check_random_matrix(Matrix, 6, 9>())); CALL_SUBTEST_13(check_random_matrix( MatrixXf(internal::random(1, EIGEN_TEST_MAX_SIZE), internal::random(1, EIGEN_TEST_MAX_SIZE)))); CALL_SUBTEST_14(check_random_matrix( MatrixXd(internal::random(1, EIGEN_TEST_MAX_SIZE), internal::random(1, EIGEN_TEST_MAX_SIZE)))); CALL_SUBTEST_15(check_random_matrix( MatrixXcf(internal::random(1, EIGEN_TEST_MAX_SIZE), internal::random(1, EIGEN_TEST_MAX_SIZE)))); CALL_SUBTEST_16(check_random_matrix( MatrixXcd(internal::random(1, EIGEN_TEST_MAX_SIZE), internal::random(1, EIGEN_TEST_MAX_SIZE)))); } }