QMCPACK
benchmark_cuSolverInverter.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) 2022 QMCPACK developers.
6 //
7 // File developed by: Mark Dewing, mdewing@anl.gov, Argonne National Laboratory
8 //
9 // File created by: Mark Dewing, mdewing@anl.gov, Argonne National Laboratory
10 //////////////////////////////////////////////////////////////////////////////////////
11 
12 #include "catch.hpp"
13 
14 #include "Configuration.h"
20 
21 namespace qmcplusplus
22 {
23 
24 TEST_CASE("cuSolverInverter_bench", "[wavefunction][benchmark]")
25 {
26 #ifdef QMC_COMPLEX
27  using FullPrecValue = std::complex<double>;
28 #else
29  using FullPrecValue = double;
30 #endif
31 
33 
34  const int N = 1024;
35 
37  Matrix<FullPrecValue> m_invT(N, N);
38  Matrix<FullPrecValue> m_invT_CPU(N, N);
40  std::complex<double> log_value;
41  m.resize(N, N);
42  m_invT.resize(N, N);
43  m_invT_CPU.resize(N, N);
44  m_invGPU.resize(N, N);
45 
48 
49  BENCHMARK("cuSolverInverter") { solver.invert_transpose(m, m_invT, m_invGPU, log_value); };
50 
52  BENCHMARK("CPU") { dmat.invert_transpose(m, m_invT_CPU, log_value); };
53 
54  auto check_matrix_result = checkMatrix(m_invT, m_invT_CPU);
55  CHECKED_ELSE(check_matrix_result.result) { FAIL(check_matrix_result.result_message); }
56 }
57 
58 } // namespace qmcplusplus
helper functions for EinsplineSetBuilder
Definition: Configuration.h:43
TEST_CASE("complex_helper", "[type_traits]")
CHECKED_ELSE(check_matrix_result.result)
std::enable_if_t< std::is_same< TMAT, T_FP >::value > invert_transpose(const Matrix< TMAT > &logdetT, Matrix< TMAT > &Ainv, Matrix< TMAT, CUDAAllocator< TMAT >> &Ainv_gpu, std::complex< TREAL > &log_value)
compute the inverse of the transpose of matrix A and its determinant value in log when T_FP and TMAT ...
auto check_matrix_result
void resize(size_type n, size_type m)
Resize the container.
Definition: OhmmsMatrix.h:99
void makeRngSpdMatrix(testing::RandomForTest< RngValueType< T >> &rng, Matrix< T > &mat_spd)
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
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
implements matrix inversion via cuSolverDN
Functor to provide scope for rng when making SpdMatrix for testing.