blob: 873845f7a0033c94e575a760b6797ff8ef7c7cba [file] [log] [blame]
Kolja Brix58e086b2021-08-23 16:00:05 +00001// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2021 Kolja Brix <kolja.brix@rwth-aachen.de>
5//
6// 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/.
9
10#include "main.h"
11#include <Eigen/SVD>
12
13
14template<typename MatrixType>
15void check_generateRandomUnitaryMatrix(const Index dim)
16{
17 const MatrixType Q = generateRandomUnitaryMatrix<MatrixType>(dim);
18
19 // validate dimensions
20 VERIFY_IS_EQUAL(Q.rows(), dim);
21 VERIFY_IS_EQUAL(Q.cols(), dim);
22
23 VERIFY_IS_UNITARY(Q);
24}
25
26template<typename VectorType, typename RealScalarType>
27void check_setupRandomSvs(const Index dim, const RealScalarType max)
28{
29 const VectorType v = setupRandomSvs<VectorType, RealScalarType>(dim, max);
30
31 // validate dimensions
32 VERIFY_IS_EQUAL(v.size(), dim);
33
34 // check entries
35 for(Index i = 0; i < v.size(); ++i)
36 VERIFY_GE(v(i), 0);
37 for(Index i = 0; i < v.size()-1; ++i)
38 VERIFY_GE(v(i), v(i+1));
39}
40
41template<typename VectorType, typename RealScalarType>
42void check_setupRangeSvs(const Index dim, const RealScalarType min, const RealScalarType max)
43{
44 const VectorType v = setupRangeSvs<VectorType, RealScalarType>(dim, min, max);
45
46 // validate dimensions
47 VERIFY_IS_EQUAL(v.size(), dim);
48
49 // check entries
50 if(dim == 1) {
51 VERIFY_IS_APPROX(v(0), min);
52 } else {
53 VERIFY_IS_APPROX(v(0), max);
54 VERIFY_IS_APPROX(v(dim-1), min);
55 }
56 for(Index i = 0; i < v.size()-1; ++i)
57 VERIFY_GE(v(i), v(i+1));
58}
59
60template<typename MatrixType, typename RealScalar, typename RealVectorType>
61void check_generateRandomMatrixSvs(const Index rows, const Index cols, const Index diag_size,
62 const RealScalar min_svs, const RealScalar max_svs)
63{
64 RealVectorType svs = setupRangeSvs<RealVectorType, RealScalar>(diag_size, min_svs, max_svs);
65
66 MatrixType M;
67 generateRandomMatrixSvs(svs, rows, cols, M);
68
69 // validate dimensions
70 VERIFY_IS_EQUAL(M.rows(), rows);
71 VERIFY_IS_EQUAL(M.cols(), cols);
72 VERIFY_IS_EQUAL(svs.size(), diag_size);
73
74 // validate singular values
75 Eigen::JacobiSVD<MatrixType> SVD(M);
76 VERIFY_IS_APPROX(svs, SVD.singularValues());
77}
78
79template<typename MatrixType>
80void check_random_matrix(const MatrixType &m)
81{
82 enum {
83 Rows = MatrixType::RowsAtCompileTime,
84 Cols = MatrixType::ColsAtCompileTime,
Erik Schultheisc20e9082021-12-10 19:27:01 +000085 DiagSize = internal::min_size_prefer_dynamic(Rows, Cols)
Kolja Brix58e086b2021-08-23 16:00:05 +000086 };
87 typedef typename MatrixType::Scalar Scalar;
88 typedef typename NumTraits<Scalar>::Real RealScalar;
89 typedef Matrix<RealScalar, DiagSize, 1> RealVectorType;
90
91 const Index rows = m.rows(), cols = m.cols();
92 const Index diag_size = (std::min)(rows, cols);
93 const RealScalar min_svs = 1.0, max_svs = 1000.0;
94
95 // check generation of unitary random matrices
96 typedef Matrix<Scalar, Rows, Rows> MatrixAType;
97 typedef Matrix<Scalar, Cols, Cols> MatrixBType;
98 check_generateRandomUnitaryMatrix<MatrixAType>(rows);
99 check_generateRandomUnitaryMatrix<MatrixBType>(cols);
100
101 // test generators for singular values
102 check_setupRandomSvs<RealVectorType, RealScalar>(diag_size, max_svs);
103 check_setupRangeSvs<RealVectorType, RealScalar>(diag_size, min_svs, max_svs);
104
105 // check generation of random matrices
106 check_generateRandomMatrixSvs<MatrixType, RealScalar, RealVectorType>(rows, cols, diag_size, min_svs, max_svs);
107}
108
109EIGEN_DECLARE_TEST(random_matrix)
110{
111 for(int i = 0; i < g_repeat; i++) {
112 CALL_SUBTEST_1(check_random_matrix(Matrix<float, 1, 1>()));
113 CALL_SUBTEST_2(check_random_matrix(Matrix<float, 4, 4>()));
114 CALL_SUBTEST_3(check_random_matrix(Matrix<float, 2, 3>()));
115 CALL_SUBTEST_4(check_random_matrix(Matrix<float, 7, 4>()));
116
117 CALL_SUBTEST_5(check_random_matrix(Matrix<double, 1, 1>()));
118 CALL_SUBTEST_6(check_random_matrix(Matrix<double, 6, 6>()));
119 CALL_SUBTEST_7(check_random_matrix(Matrix<double, 5, 3>()));
120 CALL_SUBTEST_8(check_random_matrix(Matrix<double, 4, 9>()));
121
122 CALL_SUBTEST_9(check_random_matrix(Matrix<std::complex<float>, 12, 12>()));
123 CALL_SUBTEST_10(check_random_matrix(Matrix<std::complex<float>, 7, 14>()));
124 CALL_SUBTEST_11(check_random_matrix(Matrix<std::complex<double>, 15, 11>()));
125 CALL_SUBTEST_12(check_random_matrix(Matrix<std::complex<double>, 6, 9>()));
126
127 CALL_SUBTEST_13(check_random_matrix(
128 MatrixXf(internal::random<int>(1, EIGEN_TEST_MAX_SIZE), internal::random<int>(1, EIGEN_TEST_MAX_SIZE))));
129 CALL_SUBTEST_14(check_random_matrix(
130 MatrixXd(internal::random<int>(1, EIGEN_TEST_MAX_SIZE), internal::random<int>(1, EIGEN_TEST_MAX_SIZE))));
131 CALL_SUBTEST_15(check_random_matrix(
132 MatrixXcf(internal::random<int>(1, EIGEN_TEST_MAX_SIZE), internal::random<int>(1, EIGEN_TEST_MAX_SIZE))));
133 CALL_SUBTEST_16(check_random_matrix(
134 MatrixXcd(internal::random<int>(1, EIGEN_TEST_MAX_SIZE), internal::random<int>(1, EIGEN_TEST_MAX_SIZE))));
135 }
136}