13#ifndef TLAPACK_TESTUTILS_HH
14#define TLAPACK_TESTUTILS_HH
34#ifndef TLAPACK_BUILD_STANDALONE_TESTS
35 #include <catch2/catch_template_test_macros.hpp>
36 #include <catch2/generators/catch_generators.hpp>
39 #define SKIP_TEST return
45 #define SKIP_TEST return 0
48 #define GET_FIRST_ARG(arg1, ...) arg1
52 #define DEPAREN(X) ESC(ISH X)
53 #define ISH(...) ISH __VA_ARGS__
54 #define ESC(...) ESC_(__VA_ARGS__)
55 #define ESC_(...) VAN##__VA_ARGS__
61 std::string return_scanf(
const char*);
62 std::string return_scanf(std::string);
67 if constexpr (std::is_integral<T>::value) {
72 else if constexpr (std::is_enum<T>::value) {
79 std::cout <<
"Include more cases here!\n";
82 template <
class T,
class U>
83 pair<T, U> return_scanf(pair<T, U>)
85 const T x = return_scanf(T());
86 const U y = return_scanf(U());
87 return pair<T, U>(x, y);
89 template <
class... Ts>
90 std::tuple<Ts...> return_scanf(std::tuple<Ts...>)
93 constexpr size_t N = std::tuple_size<std::tuple<Ts...>>::value;
95 std::get<0>(t) = return_scanf(std::get<0>(std::tuple<Ts...>()));
97 std::get<1>(t) = return_scanf(std::get<1>(std::tuple<Ts...>()));
99 std::get<2>(t) = return_scanf(std::get<2>(std::tuple<Ts...>()));
101 std::get<3>(t) = return_scanf(std::get<3>(std::tuple<Ts...>()));
103 std::get<4>(t) = return_scanf(std::get<4>(std::tuple<Ts...>()));
105 std::get<5>(t) = return_scanf(std::get<5>(std::tuple<Ts...>()));
106 if constexpr (N > 6) {
108 std::cout <<
"Include more cases here!\n";
115 #define TEMPLATE_TEST_CASE(TITLE, TAGS, ...) \
116 using TestType = DEPAREN(GET_FIRST_ARG(__VA_ARGS__)); \
117 int main(const int argc, const char* argv[])
119 #define GENERATE(...) \
120 tlapack::catch2::return_scanf(GET_FIRST_ARG(__VA_ARGS__))
122 #define DYNAMIC_SECTION(...) std::cout << __VA_ARGS__ << std::endl;
124 #define REQUIRE(cond) \
125 std::cout << #cond << ": " \
126 << (static_cast<bool>(cond) ? "true" : "false") << std::endl
128 #define CHECK(cond) \
129 std::cout << #cond << ": " \
130 << (static_cast<bool>(cond) ? "true" : "false") << std::endl
132 #define INFO(...) std::cout << __VA_ARGS__ << std::endl;
133 #define UNSCOPED_INFO(...) std::cout << __VA_ARGS__ << std::endl;
148template <TLAPACK_MATRIX matrix_t>
155 const idx_t m = nrows(
Q);
156 const idx_t n = ncols(
Q);
184template <TLAPACK_MATRIX matrix_t>
193 const idx_t m = min(nrows(
Q), ncols(
Q));
212template <TLAPACK_MATRIX matrix_t>
248template <TLAPACK_MATRIX matrix_t>
259 const idx_t n = ncols(
A);
263 std::vector<T>
work_;
282template <TLAPACK_MATRIX matrix_t>
326template <TLAPACK_MATRIX matrix_t>
338 const idx_t n = ncols(
A);
342 std::vector<T>
work_;
355 const LegacyMatrix<std::complex<float>,
size_t, Layout::ColMajor>& A);
357 const LegacyMatrix<std::complex<double>,
size_t, Layout::ColMajor>& A);
358void print_rowmajormatrix_r(
360void print_rowmajormatrix_d(
362void print_rowmajormatrix_c(
363 const LegacyMatrix<std::complex<float>,
size_t, Layout::RowMajor>& A);
364void print_rowmajormatrix_z(
365 const LegacyMatrix<std::complex<double>,
size_t, Layout::RowMajor>& A);
371std::string visualize_matrix_r(
373std::string visualize_matrix_d(
375std::string visualize_matrix_c(
376 const LegacyMatrix<std::complex<float>,
size_t, Layout::ColMajor>& A);
377std::string visualize_matrix_z(
378 const LegacyMatrix<std::complex<double>,
size_t, Layout::ColMajor>& A);
379std::string visualize_rowmajormatrix_r(
381std::string visualize_rowmajormatrix_d(
383std::string visualize_rowmajormatrix_c(
384 const LegacyMatrix<std::complex<float>,
size_t, Layout::RowMajor>& A);
385std::string visualize_rowmajormatrix_z(
386 const LegacyMatrix<std::complex<double>,
size_t, Layout::RowMajor>& A);
MaxtrixMarket class and random generators.
constexpr internal::FrobNorm FROB_NORM
Frobenius norm of matrices.
Definition types.hpp:347
constexpr internal::UpperTriangle UPPER_TRIANGLE
Upper Triangle access.
Definition types.hpp:186
constexpr internal::GeneralAccess GENERAL
General access.
Definition types.hpp:180
constexpr internal::ConjTranspose CONJ_TRANS
conjugate transpose
Definition types.hpp:264
constexpr internal::NoTranspose NO_TRANS
no transpose
Definition types.hpp:260
auto lange(norm_t normType, const matrix_t &A)
Calculates the norm of a matrix.
Definition lange.hpp:38
void laset(uplo_t uplo, const type_t< matrix_t > &alpha, const type_t< matrix_t > &beta, matrix_t &A)
Initializes a matrix to diagonal and off-diagonal values.
Definition laset.hpp:38
real_type< type_t< matrix_t > > check_generalized_similarity_transform(matrix_t &A, matrix_t &Q, matrix_t &Z, matrix_t &B, matrix_t &res, matrix_t &work)
Calculates res = Q'*A*Z - B and the frobenius norm of res.
Definition testutils.hpp:283
real_type< type_t< matrix_t > > check_orthogonality(matrix_t &Q, matrix_t &res)
Calculates res = Q'*Q - I if m <= n or res = Q*Q' otherwise Also computes the frobenius norm of res.
Definition testutils.hpp:149
void lacpy(uplo_t uplo, const matrixA_t &A, matrixB_t &B)
Copies a matrix from A to B.
Definition lacpy.hpp:38
auto lanhe(norm_t normType, uplo_t uplo, const matrix_t &A)
Calculates the norm of a hermitian matrix.
Definition lanhe.hpp:43
real_type< type_t< matrix_t > > check_similarity_transform(matrix_t &A, matrix_t &Q, matrix_t &B, matrix_t &res, matrix_t &work)
Calculates res = Q'*A*Q - B and the frobenius norm of res.
Definition testutils.hpp:213
void herk(Uplo uplo, Op trans, const alpha_t &alpha, const matrixA_t &A, const beta_t &beta, matrixC_t &C)
Hermitian rank-k update:
Definition herk.hpp:68
void gemm(Op transA, Op transB, const alpha_t &alpha, const matrixA_t &A, const matrixB_t &B, const beta_t &beta, matrixC_t &C)
General matrix-matrix multiply:
Definition gemm.hpp:61
#define tlapack_check(cond)
Throw an error if cond is false.
Definition exceptionHandling.hpp:98
Concept for matrices that can be converted to a legacy matrix.
typename traits::real_type_traits< Types..., int >::type real_type
The common real type of the list of types.
Definition scalar_type_traits.hpp:113
Definitions for the unit tests.