QMCPACK
qmcplusplus::cuBLAS Namespace Reference

interface to cuBLAS calls for different data types S/C/D/Z More...

Functions

cublasOperation_t convertOperation (const char trans)
 
cublasStatus_t geam (cublasHandle_t &handle, cublasOperation_t &transa, cublasOperation_t &transb, int m, int n, const float *alpha, const float *A, int lda, const float *beta, const float *B, int ldb, float *C, int ldc)
 
cublasStatus_t geam (cublasHandle_t &handle, cublasOperation_t transa, cublasOperation_t transb, int m, int n, const double *alpha, const double *A, int lda, const double *beta, const double *B, int ldb, double *C, int ldc)
 
cublasStatus_t geam (cublasHandle_t &handle, cublasOperation_t transa, cublasOperation_t transb, int m, int n, const std::complex< double > *alpha, const std::complex< double > *A, int lda, const std::complex< double > *beta, const std::complex< double > *B, int ldb, std::complex< double > *C, int ldc)
 
cublasStatus_t geam (cublasHandle_t &handle, cublasOperation_t transa, cublasOperation_t transb, int m, int n, const std::complex< float > *alpha, const std::complex< float > *A, int lda, const std::complex< float > *beta, const std::complex< float > *B, int ldb, std::complex< float > *C, int ldc)
 
cublasStatus_t getrf_batched (cublasHandle_t &handle, int n, float *A[], int lda, int *PivotArray, int *infoArray, int batchSize)
 
cublasStatus_t getrf_batched (cublasHandle_t &handle, int n, double *A[], int lda, int *PivotArray, int *infoArray, int batchSize)
 
cublasStatus_t getrf_batched (cublasHandle_t &handle, int n, std::complex< float > *A[], int lda, int *PivotArray, int *infoArray, int batchSize)
 
cublasStatus_t getrf_batched (cublasHandle_t &handle, int n, std::complex< double > *A[], int lda, int *PivotArray, int *infoArray, int batchSize)
 
cublasStatus_t getri_batched (cublasHandle_t &handle, int n, float *A[], int lda, int *PivotArray, float *C[], int ldc, int *infoArray, int batchSize)
 
cublasStatus_t getri_batched (cublasHandle_t &handle, int n, double *A[], int lda, int *PivotArray, double *C[], int ldc, int *infoArray, int batchSize)
 
cublasStatus_t getri_batched (cublasHandle_t &handle, int n, std::complex< float > *A[], int lda, int *PivotArray, std::complex< float > *C[], int ldc, int *infoArray, int batchSize)
 
cublasStatus_t getri_batched (cublasHandle_t &handle, int n, std::complex< double > *A[], int lda, int *PivotArray, std::complex< double > *C[], int ldc, int *infoArray, int batchSize)
 

Detailed Description

interface to cuBLAS calls for different data types S/C/D/Z

Function Documentation

◆ convertOperation()

cublasOperation_t qmcplusplus::cuBLAS::convertOperation ( const char  trans)
inline

Definition at line 96 of file cuBLAS.hpp.

References CUBLAS_OP_C, CUBLAS_OP_N, and CUBLAS_OP_T.

Referenced by qmcplusplus::compute::BLAS::gemm(), qmcplusplus::compute::BLAS::gemm_batched(), and qmcplusplus::compute::BLAS::gemv().

97 {
98  if (trans == 'N' || trans == 'n')
99  return CUBLAS_OP_N;
100  else if (trans == 'T' || trans == 't')
101  return CUBLAS_OP_T;
102  else if (trans == 'C' || trans == 'c')
103  return CUBLAS_OP_C;
104  else
105  throw std::runtime_error(
106  "cuBLAS::convertOperation trans can only be 'N', 'T', 'C', 'n', 't', 'c'. Input value is " +
107  std::string(1, trans));
108 }
#define CUBLAS_OP_N
Definition: cuda2hip.h:19
#define CUBLAS_OP_C
Definition: cuda2hip.h:21
#define CUBLAS_OP_T
Definition: cuda2hip.h:20

◆ geam() [1/4]

cublasStatus_t qmcplusplus::cuBLAS::geam ( cublasHandle_t handle,
cublasOperation_t transa,
cublasOperation_t transb,
int  m,
int  n,
const float *  alpha,
const float *  A,
int  lda,
const float *  beta,
const float *  B,
int  ldb,
float *  C,
int  ldc 
)
inline

Definition at line 110 of file cuBLAS.hpp.

References qmcplusplus::Units::distance::A, B(), qmcplusplus::Units::charge::C, cublasSgeam, qmcplusplus::lda, qmcplusplus::Units::distance::m, and qmcplusplus::n.

Referenced by DiracMatrixComputeCUDA< VALUE_FP >::mw_computeInvertAndLog(), and qmcplusplus::TEST_CASE().

123 {
124  return cublasSgeam(handle, transa, transb, m, n, alpha, A, lda, beta, B, ldb, C, ldc);
125 }
#define cublasSgeam
Definition: cuda2hip.h:56
double B(double x, int k, int i, const std::vector< double > &t)

◆ geam() [2/4]

cublasStatus_t qmcplusplus::cuBLAS::geam ( cublasHandle_t handle,
cublasOperation_t  transa,
cublasOperation_t  transb,
int  m,
int  n,
const double *  alpha,
const double *  A,
int  lda,
const double *  beta,
const double *  B,
int  ldb,
double *  C,
int  ldc 
)
inline

Definition at line 127 of file cuBLAS.hpp.

References qmcplusplus::Units::distance::A, B(), qmcplusplus::Units::charge::C, cublasDgeam, qmcplusplus::lda, qmcplusplus::Units::distance::m, and qmcplusplus::n.

140 {
141  return cublasDgeam(handle, transa, transb, m, n, alpha, A, lda, beta, B, ldb, C, ldc);
142 }
#define cublasDgeam
Definition: cuda2hip.h:49
double B(double x, int k, int i, const std::vector< double > &t)

◆ geam() [3/4]

cublasStatus_t qmcplusplus::cuBLAS::geam ( cublasHandle_t handle,
cublasOperation_t  transa,
cublasOperation_t  transb,
int  m,
int  n,
const std::complex< double > *  alpha,
const std::complex< double > *  A,
int  lda,
const std::complex< double > *  beta,
const std::complex< double > *  B,
int  ldb,
std::complex< double > *  C,
int  ldc 
)
inline

Definition at line 144 of file cuBLAS.hpp.

References qmcplusplus::Units::distance::A, B(), qmcplusplus::Units::charge::C, castNativeType, cublasZgeam, qmcplusplus::lda, qmcplusplus::Units::distance::m, and qmcplusplus::n.

157 {
158  return cublasZgeam(handle, transa, transb, m, n, castNativeType(alpha), castNativeType(A), lda, castNativeType(beta),
159  castNativeType(B), ldb, castNativeType(C), ldc);
160 }
#define castNativeType
Definition: cuBLAS.hpp:23
#define cublasZgeam
Definition: cuda2hip.h:63
double B(double x, int k, int i, const std::vector< double > &t)

◆ geam() [4/4]

cublasStatus_t qmcplusplus::cuBLAS::geam ( cublasHandle_t handle,
cublasOperation_t  transa,
cublasOperation_t  transb,
int  m,
int  n,
const std::complex< float > *  alpha,
const std::complex< float > *  A,
int  lda,
const std::complex< float > *  beta,
const std::complex< float > *  B,
int  ldb,
std::complex< float > *  C,
int  ldc 
)
inline

Definition at line 162 of file cuBLAS.hpp.

References qmcplusplus::Units::distance::A, B(), qmcplusplus::Units::charge::C, castNativeType, cublasCgeam, qmcplusplus::lda, qmcplusplus::Units::distance::m, and qmcplusplus::n.

175 {
176  return cublasCgeam(handle, transa, transb, m, n, castNativeType(alpha), castNativeType(A), lda, castNativeType(beta),
177  castNativeType(B), ldb, castNativeType(C), ldc);
178 }
#define cublasCgeam
Definition: cuda2hip.h:42
#define castNativeType
Definition: cuBLAS.hpp:23
double B(double x, int k, int i, const std::vector< double > &t)

◆ getrf_batched() [1/4]

cublasStatus_t qmcplusplus::cuBLAS::getrf_batched ( cublasHandle_t handle,
int  n,
float *  A[],
int  lda,
int *  PivotArray,
int *  infoArray,
int  batchSize 
)
inline

Definition at line 180 of file cuBLAS.hpp.

References qmcplusplus::Units::distance::A, cublasSgetrfBatched, qmcplusplus::lda, and qmcplusplus::n.

187 {
188  return cublasSgetrfBatched(handle, n, A, lda, PivotArray, infoArray, batchSize);
189 }
#define cublasSgetrfBatched
Definition: cuda2hip.h:61

◆ getrf_batched() [2/4]

cublasStatus_t qmcplusplus::cuBLAS::getrf_batched ( cublasHandle_t handle,
int  n,
double *  A[],
int  lda,
int *  PivotArray,
int *  infoArray,
int  batchSize 
)
inline

Definition at line 191 of file cuBLAS.hpp.

References qmcplusplus::Units::distance::A, cublasDgetrfBatched, qmcplusplus::lda, and qmcplusplus::n.

198 {
199  return cublasDgetrfBatched(handle, n, A, lda, PivotArray, infoArray, batchSize);
200 }
#define cublasDgetrfBatched
Definition: cuda2hip.h:54

◆ getrf_batched() [3/4]

cublasStatus_t qmcplusplus::cuBLAS::getrf_batched ( cublasHandle_t handle,
int  n,
std::complex< float > *  A[],
int  lda,
int *  PivotArray,
int *  infoArray,
int  batchSize 
)
inline

Definition at line 202 of file cuBLAS.hpp.

References qmcplusplus::Units::distance::A, qmcplusplus::castCUDAType(), cublasCgetrfBatched, qmcplusplus::lda, and qmcplusplus::n.

209 {
210  return cublasCgetrfBatched(handle, n, castCUDAType(A), lda, PivotArray, infoArray, batchSize);
211 }
#define cublasCgetrfBatched
Definition: cuda2hip.h:47
CUDATypeMap< T > castCUDAType(T var)

◆ getrf_batched() [4/4]

cublasStatus_t qmcplusplus::cuBLAS::getrf_batched ( cublasHandle_t handle,
int  n,
std::complex< double > *  A[],
int  lda,
int *  PivotArray,
int *  infoArray,
int  batchSize 
)
inline

Definition at line 213 of file cuBLAS.hpp.

References qmcplusplus::Units::distance::A, qmcplusplus::castCUDAType(), cublasZgetrfBatched, qmcplusplus::lda, and qmcplusplus::n.

220 {
221  return cublasZgetrfBatched(handle, n, castCUDAType(A), lda, PivotArray, infoArray, batchSize);
222 }
CUDATypeMap< T > castCUDAType(T var)
#define cublasZgetrfBatched
Definition: cuda2hip.h:68

◆ getri_batched() [1/4]

cublasStatus_t qmcplusplus::cuBLAS::getri_batched ( cublasHandle_t handle,
int  n,
float *  A[],
int  lda,
int *  PivotArray,
float *  C[],
int  ldc,
int *  infoArray,
int  batchSize 
)
inline

Definition at line 224 of file cuBLAS.hpp.

References qmcplusplus::Units::distance::A, qmcplusplus::Units::charge::C, cublasSgetriBatched, qmcplusplus::lda, and qmcplusplus::n.

233 {
234  return cublasSgetriBatched(handle, n, A, lda, PivotArray, C, ldc, infoArray, batchSize);
235 }
#define cublasSgetriBatched
Definition: cuda2hip.h:62

◆ getri_batched() [2/4]

cublasStatus_t qmcplusplus::cuBLAS::getri_batched ( cublasHandle_t handle,
int  n,
double *  A[],
int  lda,
int *  PivotArray,
double *  C[],
int  ldc,
int *  infoArray,
int  batchSize 
)
inline

Definition at line 237 of file cuBLAS.hpp.

References qmcplusplus::Units::distance::A, qmcplusplus::Units::charge::C, cublasDgetriBatched, qmcplusplus::lda, and qmcplusplus::n.

246 {
247  return cublasDgetriBatched(handle, n, A, lda, PivotArray, C, ldc, infoArray, batchSize);
248 }
#define cublasDgetriBatched
Definition: cuda2hip.h:55

◆ getri_batched() [3/4]

cublasStatus_t qmcplusplus::cuBLAS::getri_batched ( cublasHandle_t handle,
int  n,
std::complex< float > *  A[],
int  lda,
int *  PivotArray,
std::complex< float > *  C[],
int  ldc,
int *  infoArray,
int  batchSize 
)
inline

Definition at line 250 of file cuBLAS.hpp.

References qmcplusplus::Units::distance::A, qmcplusplus::Units::charge::C, qmcplusplus::castCUDAType(), cublasCgetriBatched, qmcplusplus::lda, and qmcplusplus::n.

259 {
260  return cublasCgetriBatched(handle, n, castCUDAType(A), lda, PivotArray, castCUDAType(C), ldc, infoArray, batchSize);
261 }
#define cublasCgetriBatched
Definition: cuda2hip.h:48
CUDATypeMap< T > castCUDAType(T var)

◆ getri_batched() [4/4]

cublasStatus_t qmcplusplus::cuBLAS::getri_batched ( cublasHandle_t handle,
int  n,
std::complex< double > *  A[],
int  lda,
int *  PivotArray,
std::complex< double > *  C[],
int  ldc,
int *  infoArray,
int  batchSize 
)
inline

Definition at line 263 of file cuBLAS.hpp.

References qmcplusplus::Units::distance::A, qmcplusplus::Units::charge::C, qmcplusplus::castCUDAType(), cublasZgetriBatched, qmcplusplus::lda, and qmcplusplus::n.

272 {
273  return cublasZgetriBatched(handle, n, castCUDAType(A), lda, PivotArray, castCUDAType(C), ldc, infoArray, batchSize);
274 }
#define cublasZgetriBatched
Definition: cuda2hip.h:69
CUDATypeMap< T > castCUDAType(T var)