ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
math_kernel_op.h
Go to the documentation of this file.
1// TODO: This is a temperary location for these functions.
2// And will be moved to a global module(module base) later.
3#ifndef MODULE_HSOLVER_MATH_KERNEL_H
4#define MODULE_HSOLVER_MATH_KERNEL_H
5
7
8
11
12#if defined(__CUDA) || defined(__UT_USE_CUDA)
13#include <cuda_runtime.h>
14
15#include "cublas_v2.h"
16#endif //__CUDA || __UT_USE_CUDA
17
18namespace ModuleBase {
19
20//---------------------------------------------------------------------------------
21//-----------------------------0. Tool Functions-----------------------------------
22//---------------------------------------------------------------------------------
23inline std::complex<double> set_real_tocomplex(const std::complex<double> &x) {
24 return {x.real(), 0.0};
25}
26
27inline std::complex<float> set_real_tocomplex(const std::complex<float> &x) {
28 return {x.real(), 0.0};
29}
30
31inline double set_real_tocomplex(const double &x) { return x; }
32
33inline float set_real_tocomplex(const float &x) { return x; }
34
35inline std::complex<double> get_conj(const std::complex<double> &x) {
36 return {x.real(), -x.imag()};
37}
38
39inline std::complex<float> get_conj(const std::complex<float> &x) {
40 return {x.real(), -x.imag()};
41}
42
43inline double get_conj(const double &x) { return x; }
44
45inline float get_conj(const float &x) { return x; }
46
47
48//---------------------------------------------------------------------------------
49//-----------------------------1. Vector Operations--------------------------------
50//---------------------------------------------------------------------------------
51template <typename FPTYPE, typename Device> struct scal_op {
62 void operator()(const int &N,
63 const std::complex<FPTYPE> *alpha, std::complex<FPTYPE> *X,
64 const int &incx);
65};
66
67template <typename T, typename Device> struct vector_mul_real_op {
68 using Real = typename GetTypeReal<T>::type;
80 void operator()(const int dim, T* result, const T* vector, const Real constant);
81};
82
83// vector operator: result[i] = vector1[i](complex) * vector2[i](not complex)
84template <typename T, typename Device> struct vector_mul_vector_op {
85 using Real = typename GetTypeReal<T>::type;
96 void operator()(const int& dim, T* result, const T* vector1, const Real* vector2, const bool& add = false);
97};
98
99// vector operator: result[i] = vector[i] / constant
100template <typename T, typename Device> struct vector_div_constant_op {
101 using Real = typename GetTypeReal<T>::type;
111 void operator()(const int& dim, T* result, const T* vector, const Real constant);
112};
113
114// vector operator: result[i] = vector1[i](complex) / vector2[i](not complex)
115template <typename T, typename Device> struct vector_div_vector_op {
116 using Real = typename GetTypeReal<T>::type;
126 void operator()(const int &dim, T *result, const T *vector1,
127 const Real *vector2);
128};
129
130// compute Y = alpha * X + Y
131template <typename T, typename Device> struct axpy_op {
144 void operator()(const int &N, const T *alpha, const T *X,
145 const int &incX, T *Y, const int &incY);
146};
147
148// vector operator: result[i] = vector1[i] * constant1 + vector2[i] * constant2
149template <typename T, typename Device>
151 using Real = typename GetTypeReal<T>::type;
163 void operator()(const int &dim, T *result, const T *vector1,
164 const Real constant1, const T *vector2, const Real constant2);
165};
166
167template <typename T, typename Device> struct dot_real_op {
168 using Real = typename GetTypeReal<T>::type;
181 Real operator()(const int &dim, const T *psi_L,
182 const T *psi_R, const bool reduce = true);
183};
184
185
186//---------------------------------------------------------------------------------
187//-----------------------------2. Matrix Operations--------------------------------
188//---------------------------------------------------------------------------------
189
190// compute y = alpha * op(A) * x + beta * y
191template <typename T, typename Device> struct gemv_op {
209 void operator()(const char &trans, const int &m,
210 const int &n, const T *alpha, const T *A, const int &lda,
211 const T *X, const int &incx, const T *beta, T *Y,
212 const int &incy);
213};
214
215// compute C = alpha * op(A) * op(B) + beta * C
216template <typename T, typename Device> struct gemm_op {
236 void operator()(const char &transa, const char &transb,
237 const int &m, const int &n, const int &k, const T *alpha,
238 const T *a, const int &lda, const T *b, const int &ldb,
239 const T *beta, T *c, const int &ldc);
240};
241
242#ifdef __DSP
243// compute Y = alpha * op(A) * X + beta * Y on DSP Hardware
244template <typename T, typename Device> struct gemv_op_mt {
262 void operator()(const char &trans, const int &m,
263 const int &n, const T *alpha, const T *A, const int &lda,
264 const T *X, const int &incx, const T *beta, T *Y,
265 const int &incy);
266};
267
268// compute C = alpha * op(A) * op(B) + beta * C on DSP Hardware
269template <typename T, typename Device> struct gemm_op_mt {
289 void operator()(const char &transa, const char &transb,
290 const int &m, const int &n, const int &k, const T *alpha,
291 const T *a, const int &lda, const T *b, const int &ldb,
292 const T *beta, T *c, const int &ldc);
293};
294#endif
295
296template <typename T, typename Device> struct matrixTranspose_op {
306 void operator()(const int &row, const int &col,
307 const T *input_matrix, T *output_matrix);
308};
309
310template <typename T, typename Device> struct matrixCopy {
322 void operator()(const int& n1, const int& n2, const T* A, const int& LDA, T* B, const int& LDB);
323};
324
325template <typename T, typename Device>
327 using Real = typename GetTypeReal<T>::type;
341 void operator()(const int &m, const int &n,
342 T *a,
343 const int &lda,
344 const Real *b,
345 const Real alpha,
346 T *c,
347 const int &ldc);
348};
349
350template <typename T, typename Device>
352 using Real = typename GetTypeReal<T>::type;
353
354 void operator()(const Device *d, const int &nbase, const int &nbase_x, const int &notconv,
355 T *result, const T *vectors, const Real *eigenvalues);
356};
357
358template <typename T, typename Device>
360 using Real = typename GetTypeReal<T>::type;
361 void operator()(const Device* d,
362 const int& dim,
363 T* psi_iter,
364 const int& nbase,
365 const int& notconv,
366 const Real* precondition,
367 const Real* eigenvalues);
368};
369
370template <typename T, typename Device>
372 using Real = typename GetTypeReal<T>::type;
373 void operator()(const Device* d,
374 const int& dim,
375 T* psi_iter,
376 const int& nbase,
377 const int& notconv,
378 Real* psi_norm = nullptr);
379};
380
381template <typename T>
382struct normalize_op<T, base_device::DEVICE_GPU> {
383 using Real = typename GetTypeReal<T>::type;
384 void operator()(const base_device::DEVICE_GPU* d,
385 const int& dim,
386 T* psi_iter,
387 const int& nbase,
388 const int& notconv,
389 Real* psi_norm);
390};
391
392#if __CUDA || __UT_USE_CUDA || __ROCM || __UT_USE_ROCM
393// Partially specialize functor for base_device::GpuDevice.
394template <typename T> struct dot_real_op<T, base_device::DEVICE_GPU> {
395 using Real = typename GetTypeReal<T>::type;
396 Real operator()(const int &dim,
397 const T *psi_L, const T *psi_R, const bool reduce = true);
398};
399
400// vector operator: result[i] = vector[i] / constant
401template <typename T>
402struct vector_mul_real_op<T, base_device::DEVICE_GPU>
403{
404 using Real = typename GetTypeReal<T>::type;
405 void operator()(const int dim, T* result, const T* vector, const Real constant);
406};
407
408// vector operator: result[i] = vector1[i](complex) * vector2[i](not complex)
409template <typename T> struct vector_mul_vector_op<T, base_device::DEVICE_GPU> {
410 using Real = typename GetTypeReal<T>::type;
411 void operator()(const int& dim, T* result, const T* vector1, const Real* vector2, const bool& add = false);
412};
413
414// vector operator: result[i] = vector[i] / constant
415template <typename T> struct vector_div_constant_op<T, base_device::DEVICE_GPU> {
416 using Real = typename GetTypeReal<T>::type;
417 void operator()(const int& dim, T* result, const T* vector, const Real constant);
418};
419
420// vector operator: result[i] = vector1[i](complex) / vector2[i](not complex)
421template <typename T> struct vector_div_vector_op<T, base_device::DEVICE_GPU> {
422 using Real = typename GetTypeReal<T>::type;
423 void operator()(const int &dim, T *result,
424 const T *vector1, const Real *vector2);
425};
426
427// vector operator: result[i] = vector1[i] * constant1 + vector2[i] * constant2
428template <typename T>
429struct vector_add_vector_op<T, base_device::DEVICE_GPU> {
430 using Real = typename GetTypeReal<T>::type;
431 void operator()(const int &dim, T *result,
432 const T *vector1, const Real constant1, const T *vector2,
433 const Real constant2);
434};
435
436template <typename T> struct matrixCopy<T, base_device::DEVICE_GPU> {
437 void operator()(const int& n1,
438 const int& n2,
439 const T* A, // input
440 const int& LDA,
441 T* B, // output
442 const int& LDB);
443};
444
445template <typename T> struct matrix_mul_vector_op<T, base_device::DEVICE_GPU> {
446 using Real = typename GetTypeReal<T>::type;
447 void operator()(const int &m, const int &n,
448 T *a,
449 const int &lda,
450 const Real *b,
451 const Real alpha,
452 T *c,
453 const int &ldc);
454};
455
456void createGpuBlasHandle();
457void destoryBLAShandle();
458
459// vector operator: result[i] = -lambda[i] * vector[i]
460template <typename T> struct apply_eigenvalues_op<T, base_device::DEVICE_GPU> {
461 using Real = typename GetTypeReal<T>::type;
462
463 void operator()(const base_device::DEVICE_GPU *d, const int &nbase, const int &nbase_x, const int &notconv,
464 T *result, const T *vectors, const Real *eigenvalues);
465};
466
467template <typename T>
468struct precondition_op<T, base_device::DEVICE_GPU> {
469 using Real = typename GetTypeReal<T>::type;
470 void operator()(const base_device::DEVICE_GPU* d,
471 const int& dim,
472 T* psi_iter,
473 const int& nbase,
474 const int& notconv,
475 const Real* precondition,
476 const Real* eigenvalues);
477};
478
479#endif // __CUDA || __UT_USE_CUDA || __ROCM || __UT_USE_ROCM
480} // namespace hsolver
481
482#endif // MODULE_HSOLVER_MATH_KERNEL_H
#define N
Definition exp.cpp:24
#define T
Definition exp.cpp:237
Definition clebsch_gordan_coeff.cpp:8
std::complex< double > set_real_tocomplex(const std::complex< double > &x)
Definition math_kernel_op.h:23
std::complex< double > get_conj(const std::complex< double > &x)
Definition math_kernel_op.h:35
Definition device.cpp:21
T type
Definition macros.h:8
Definition math_kernel_op.h:351
void operator()(const Device *d, const int &nbase, const int &nbase_x, const int &notconv, T *result, const T *vectors, const Real *eigenvalues)
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:352
Definition math_kernel_op.h:131
void operator()(const int &N, const T *alpha, const T *X, const int &incX, T *Y, const int &incY)
Y = alpha * X + Y.
Definition math_kernel_op.h:167
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:168
Real operator()(const int &dim, const T *psi_L, const T *psi_R, const bool reduce=true)
dot_real_op computes the dot product of the given complex arrays(treated as float arrays)....
Definition math_kernel_op.h:216
void operator()(const char &transa, const char &transb, const int &m, const int &n, const int &k, const T *alpha, const T *a, const int &lda, const T *b, const int &ldb, const T *beta, T *c, const int &ldc)
C = alpha * op(A) * op(B) + beta * C.
Definition math_kernel_op.h:191
void operator()(const char &trans, const int &m, const int &n, const T *alpha, const T *A, const int &lda, const T *X, const int &incx, const T *beta, T *Y, const int &incy)
y = alpha * op(A) * x + beta * y
Definition math_kernel_op.h:310
void operator()(const int &n1, const int &n2, const T *A, const int &LDA, T *B, const int &LDB)
copy matrix A to B, they can have different leading dimensions
Definition math_kernel_op.h:296
void operator()(const int &row, const int &col, const T *input_matrix, T *output_matrix)
transpose the input matrix
Definition math_kernel_op.h:326
void operator()(const int &m, const int &n, T *a, const int &lda, const Real *b, const Real alpha, T *c, const int &ldc)
a * b * beta by each column
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:327
void operator()(const base_device::DEVICE_GPU *d, const int &dim, T *psi_iter, const int &nbase, const int &notconv, Real *psi_norm)
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:383
Definition math_kernel_op.h:371
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:372
void operator()(const Device *d, const int &dim, T *psi_iter, const int &nbase, const int &notconv, Real *psi_norm=nullptr)
Definition math_kernel_op.h:359
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:360
void operator()(const Device *d, const int &dim, T *psi_iter, const int &nbase, const int &notconv, const Real *precondition, const Real *eigenvalues)
Definition math_kernel_op.h:51
void operator()(const int &N, const std::complex< FPTYPE > *alpha, std::complex< FPTYPE > *X, const int &incx)
x = alpha * x, where alpha and x are complex numbers
Definition math_kernel_op.h:150
void operator()(const int &dim, T *result, const T *vector1, const Real constant1, const T *vector2, const Real constant2)
result[i] = vector1[i] * constant1 + vector2[i] * constant2
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:151
Definition math_kernel_op.h:100
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:101
void operator()(const int &dim, T *result, const T *vector, const Real constant)
result[i] = vector[i] / constant
Definition math_kernel_op.h:115
void operator()(const int &dim, T *result, const T *vector1, const Real *vector2)
result[i] = vector1[i](complex) / vector2[i](not complex)
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:116
Definition math_kernel_op.h:67
void operator()(const int dim, T *result, const T *vector, const Real constant)
result[i] = vector[i] * constant, where vector is complex number and constant is real number。 It is d...
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:68
Definition math_kernel_op.h:84
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:85
void operator()(const int &dim, T *result, const T *vector1, const Real *vector2, const bool &add=false)
result[i] = vector1[i](complex) * vector2[i](not complex)