QMCPACK
AccelMatrixUpdateSYCL.hpp
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) 2024 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 
13 #ifndef QMCPLUSPLUS_COMPUTE_MATRIX_UPDATE_SYCL_H
14 #define QMCPLUSPLUS_COMPUTE_MATRIX_UPDATE_SYCL_H
15 
16 #include <QueueAliases.hpp>
17 #include "matrix_update_helper.hpp"
18 
19 namespace qmcplusplus
20 {
21 
22 namespace compute
23 {
24 
25 template<typename T>
27  const int rowchanged,
28  const int n,
29  const T* const Ainv[],
30  const int lda,
31  T* const temp[],
32  T* const rcopy[],
33  const T* const phi_vgl_in[],
34  const size_t phi_vgl_stride,
35  T* const dphi_out[],
36  T* const d2phi_out[],
37  const int batch_count)
38 {
39  try
40  {
41  SYCL::copyAinvRow_saveGL_batched(queue.getNative(), rowchanged, n, Ainv, lda, temp, rcopy, phi_vgl_in,
42  phi_vgl_stride, dphi_out, d2phi_out, batch_count);
43  }
44  catch (sycl::exception& e)
45  {
46  throw std::runtime_error(std::string("SYCL::copyAinvRow_saveGL_batched exception: ") + e.what());
47  }
48 }
49 
50 template<typename T>
52  const int n,
53  const T* const Ainvrow[],
54  const T* const dpsiMrow[],
55  T* const grads_now,
56  const int batch_count)
57 {
58  try
59  {
60  SYCL::calcGradients_batched(queue.getNative(), n, Ainvrow, dpsiMrow, grads_now, batch_count);
61  }
62  catch (sycl::exception& e)
63  {
64  throw std::runtime_error(std::string("SYCL::calcGradients_batched exception: ") + e.what());
65  }
66 }
67 
68 template<typename T>
70  int* const delay_list[],
71  const int rowchanged,
72  const int delay_count,
73  T* const binv[],
74  const int binv_lda,
75  const T* const ratio_inv,
76  const T* const phi_vgl_in[],
77  const size_t phi_vgl_stride,
78  T* const phi_out[],
79  T* const dphi_out[],
80  T* const d2phi_out[],
81  const int norb,
82  const int n_accepted,
83  const int batch_count)
84 {
85  try
86  {
87  SYCL::add_delay_list_save_sigma_VGL_batched(queue.getNative(), delay_list, rowchanged, delay_count, binv, binv_lda,
88  ratio_inv, phi_vgl_in, phi_vgl_stride, phi_out, dphi_out, d2phi_out,
89  norb, n_accepted, batch_count);
90  }
91  catch (sycl::exception& e)
92  {
93  throw std::runtime_error(std::string("SYCL::add_delay_list_save_y_VGL_batched exception: ") + e.what());
94  }
95 }
96 
97 
98 template<typename T>
100  const int* const delay_list[],
101  const int delay_count,
102  T* const tempMat[],
103  const int lda,
104  const int batch_count)
105 {
106  try
107  {
108  SYCL::applyW_batched(queue.getNative(), delay_list, delay_count, tempMat, lda, batch_count);
109  }
110  catch (sycl::exception& e)
111  {
112  throw std::runtime_error(std::string("SYCL::applyW_batched exception: ") + e.what());
113  }
114 }
115 
116 
117 } // namespace compute
118 } // namespace qmcplusplus
119 #endif
void add_delay_list_save_sigma_VGL_batched(Queue< PlatformKind::CUDA > &queue, int *const delay_list[], const int rowchanged, const int delay_count, T *const binv[], const int binv_lda, const T *const ratio_inv, const T *const phi_vgl_in[], const size_t phi_vgl_stride, T *const phi_out[], T *const dphi_out[], T *const d2phi_out[], const int norb, const int n_accepted, const int batch_count)
void applyW_batched(Queue< PlatformKind::CUDA > &queue, const int *const delay_list[], const int delay_count, T *const tempMat[], const int lda, const int batch_count)
helper functions for EinsplineSetBuilder
Definition: Configuration.h:43
sycl::event calcGradients_batched(sycl::queue &aq, const int n, const T *const Ainvrow[], const T *const dpsiMrow[], T *const grads_now, const int batch_count, const std::vector< sycl::event > &dependencies)
sycl::event add_delay_list_save_sigma_VGL_batched(sycl::queue &aq, int *const delay_list[], const int rowchanged, const int delay_count, T *const binv[], const int binv_lda, const T *const ratio_inv, const T *const phi_vgl_in[], const int phi_vgl_stride, T *const phi_out[], T *const dphi_out[], T *const d2phi_out[], const int norb, const int n_accepted, const int batch_count, const std::vector< sycl::event > &dependencies)
sycl::event applyW_batched(sycl::queue &aq, const int *const delay_list[], const int delay_count, T *const tempMat[], const int lda, const int batch_count, const std::vector< sycl::event > &dependencies)
void calcGradients_batched(Queue< PlatformKind::CUDA > &queue, const int n, const T *const Ainvrow[], const T *const dpsiMrow[], T *const grads_now, const int batch_count)
void copyAinvRow_saveGL_batched(Queue< PlatformKind::CUDA > &queue, const int rowchanged, const int n, const T *const Ainv[], const int lda, T *const temp[], T *const rcopy[], const T *const phi_vgl_in[], const size_t phi_vgl_stride, T *const dphi_out[], T *const d2phi_out[], const int batch_count)
sycl::event copyAinvRow_saveGL_batched(sycl::queue &aq, const int rowchanged, const int n, const T *const Ainv[], const int lda, T *const temp[], T *const rcopy[], const T *const phi_vgl_in[], const int phi_vgl_stride, T *const dphi_out[], T *const d2phi_out[], const int batch_count, const std::vector< sycl::event > &dependencies)