QMCPACK
test_syclSolverInverter.cpp
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////////////
2 // This file is distributed under the University of Illinois/NCSA Open Source License.
3 // See LICENSE file in top directory for details.
4 //
5 // Copyright (c) 2021 QMCPACK developers.
6 //
7 // File developed by: Ye Luo, yeluo@anl.gov, Argonne National Laboratory
8 //
9 // File created by: Ye Luo, yeluo@anl.gov, Argonne National Laboratory
10 //////////////////////////////////////////////////////////////////////////////////////
11 
12 #include "catch.hpp"
13 
14 #include "Configuration.h"
20 
21 namespace qmcplusplus
22 {
23 template<typename T, typename T_FP>
24 void test_inverse(const std::int64_t N)
25 {
27 
29 
30  Matrix<T> m(N, N);
31  Matrix<T> m_invT(N, N);
32  Matrix<T> m_invT_CPU(N, N);
34  m_invGPU.resize(N, N);
35 
36  std::complex<double> log_value, log_value_cpu;
39 
40  solver.invert_transpose(m, m_invT, m_invGPU, log_value, m_queue);
41  m_queue.wait();
42 
43  DiracMatrix<T_FP> dmat;
44  dmat.invert_transpose(m, m_invT_CPU, log_value_cpu);
45 
46  CHECK(log_value == ComplexApprox(log_value_cpu));
47 
48  auto check_matrix_result = checkMatrix(m_invT, m_invT_CPU);
49  CHECKED_ELSE(check_matrix_result.result) { FAIL(check_matrix_result.result_message); }
50 }
51 
52 TEMPLATE_TEST_CASE("syclSolverInverter", "[wavefunction][fermion]", double, float)
53 {
54  // TestType is defined by Catch. It is the type in each instantiation of the templated test case.
55 #ifdef QMC_COMPLEX
56  using FullPrecValue = std::complex<double>;
57  using Value = std::complex<TestType>;
58 #else
59  using FullPrecValue = double;
60  using Value = TestType;
61 #endif
62 
63  SECTION("N=117") { test_inverse<Value, FullPrecValue>(117); }
64 
65  SECTION("N=911") { test_inverse<Value, FullPrecValue>(911); }
66 }
67 
68 } // namespace qmcplusplus
sycl::queue & getSYCLDefaultDeviceDefaultQueue()
return a reference to the per-device default queue
Definition: SYCLruntime.cpp:18
helper functions for EinsplineSetBuilder
Definition: Configuration.h:43
CHECKED_ELSE(check_matrix_result.result)
auto check_matrix_result
void resize(size_type n, size_type m)
Resize the container.
Definition: OhmmsMatrix.h:99
implements matrix inversion via cuSolverDN
void makeRngSpdMatrix(testing::RandomForTest< RngValueType< T >> &rng, Matrix< T > &mat_spd)
helper class to compute matrix inversion and the log value of determinant
Definition: DiracMatrix.h:111
void test_inverse(const std::int64_t N)
CheckMatrixResult checkMatrix(M1 &a_mat, M2 &b_mat, const bool check_all=false, std::optional< const double > eps=std::nullopt)
This function checks equality a_mat and b_mat elements M1, M2 need to have their element type declare...
Definition: checkMatrix.hpp:63
CHECK(log_values[0]==ComplexApprox(std::complex< double >{ 5.603777579195571, -6.1586603331188225 }))
TEMPLATE_TEST_CASE("RandomRotationMatrix", "[numerics]", float, double)
std::enable_if_t< std::is_same< T_FP, TMAT >::value > invert_transpose(const Matrix< TMAT, ALLOC1 > &amat, Matrix< TMAT, ALLOC2 > &invMat, std::complex< TREAL > &LogDet)
compute the inverse of the transpose of matrix A and its determinant value in log when T_FP and TMAT ...
Definition: DiracMatrix.h:188
std::enable_if_t< std::is_same< TMAT, T_FP >::value > invert_transpose(const Matrix< TMAT > &logdetT, Matrix< TMAT > &Ainv, Matrix< TMAT, SYCLAllocator< TMAT >> &Ainv_gpu, std::complex< TREAL > &log_value, sycl::queue &m_queue)
compute the inverse of the transpose of matrix A and its determinant value in log when T_FP and TMAT ...
Functor to provide scope for rng when making SpdMatrix for testing.