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
9
12
13#if defined(__CUDA) || defined(__UT_USE_CUDA)
14#include <cuda_runtime.h>
15
16#include "cublas_v2.h"
17#endif //__CUDA || __UT_USE_CUDA
18
19namespace ModuleBase {
20
21//---------------------------------------------------------------------------------
22//-----------------------------0. Tool Functions-----------------------------------
23//---------------------------------------------------------------------------------
24inline std::complex<double> set_real_tocomplex(const std::complex<double> &x) {
25 return {x.real(), 0.0};
26}
27
28inline std::complex<float> set_real_tocomplex(const std::complex<float> &x) {
29 return {x.real(), 0.0};
30}
31
32inline double set_real_tocomplex(const double &x) { return x; }
33
34inline float set_real_tocomplex(const float &x) { return x; }
35
36inline std::complex<double> get_conj(const std::complex<double> &x) {
37 return {x.real(), -x.imag()};
38}
39
40inline std::complex<float> get_conj(const std::complex<float> &x) {
41 return {x.real(), -x.imag()};
42}
43
44inline double get_conj(const double &x) { return x; }
45
46inline float get_conj(const float &x) { return x; }
47
48
49//---------------------------------------------------------------------------------
50//-----------------------------1. Vector Operations--------------------------------
51//---------------------------------------------------------------------------------
52template <typename FPTYPE, typename Device> struct scal_op {
63 void operator()(const int &N,
64 const std::complex<FPTYPE> *alpha, std::complex<FPTYPE> *X,
65 const int &incx);
66};
67
68template <typename T, typename Device> struct vector_mul_real_op {
69 using Real = typename GetTypeReal<T>::type;
81 void operator()(const int dim, T* result, const T* vector, const Real constant);
82};
83
84// vector operator: result[i] = vector1[i](complex) * vector2[i](not complex)
85template <typename T, typename Device> struct vector_mul_vector_op {
86 using Real = typename GetTypeReal<T>::type;
97 void operator()(const int& dim, T* result, const T* vector1, const Real* vector2, const bool& add = false);
98};
99
100// vector operator: result[i] = vector[i] / constant
101template <typename T, typename Device> struct vector_div_constant_op {
102 using Real = typename GetTypeReal<T>::type;
112 void operator()(const int& dim, T* result, const T* vector, const Real constant);
113};
114
115// vector operator: result[i] = vector1[i](complex) / vector2[i](not complex)
116template <typename T, typename Device> struct vector_div_vector_op {
117 using Real = typename GetTypeReal<T>::type;
127 void operator()(const int &dim, T *result, const T *vector1,
128 const Real *vector2);
129};
130
131// compute Y = alpha * X + Y
132template <typename T, typename Device> struct axpy_op {
145 void operator()(const int &N, const T *alpha, const T *X,
146 const int &incX, T *Y, const int &incY);
147};
148
149// vector operator: result[i] = vector1[i] * constant1 + vector2[i] * constant2
150template <typename T, typename Device>
152 using Real = typename GetTypeReal<T>::type;
164 void operator()(const int &dim, T *result, const T *vector1,
165 const Real constant1, const T *vector2, const Real constant2);
166};
167
168template <typename T, typename Device> struct dot_real_op {
169 using Real = typename GetTypeReal<T>::type;
182 Real operator()(const int &dim, const T *psi_L,
183 const T *psi_R, const bool reduce = true);
184};
185
186
187//---------------------------------------------------------------------------------
188//-----------------------------2. Matrix Operations--------------------------------
189//---------------------------------------------------------------------------------
190
191// compute y = alpha * op(A) * x + beta * y
192template <typename T, typename Device> struct gemv_op {
210 void operator()(const char &trans, const int &m,
211 const int &n, const T *alpha, const T *A, const int &lda,
212 const T *X, const int &incx, const T *beta, T *Y,
213 const int &incy);
214};
215
216// compute C = alpha * op(A) * op(B) + beta * C
217template <typename T, typename Device> struct gemm_op {
237 void operator()(const char &transa, const char &transb,
238 const int &m, const int &n, const int &k, const T *alpha,
239 const T *a, const int &lda, const T *b, const int &ldb,
240 const T *beta, T *c, const int &ldc);
241};
242
243#ifdef __DSP
244// compute C = alpha * op(A) * op(B) + beta * C on DSP Hardware
245template <typename T, typename Device> struct gemm_op_mt {
265 void operator()(const char &transa, const char &transb,
266 const int &m, const int &n, const int &k, const T *alpha,
267 const T *a, const int &lda, const T *b, const int &ldb,
268 const T *beta, T *c, const int &ldc);
269};
270#endif
271
272template <typename T, typename Device> struct matrixTranspose_op {
282 void operator()(const int &row, const int &col,
283 const T *input_matrix, T *output_matrix);
284};
285
286template <typename T, typename Device> struct matrixCopy {
298 void operator()(const int& n1, const int& n2, const T* A, const int& LDA, T* B, const int& LDB);
299};
300
301template <typename T, typename Device>
303 using Real = typename GetTypeReal<T>::type;
317 void operator()(const int &m, const int &n,
318 T *a,
319 const int &lda,
320 const Real *b,
321 const Real alpha,
322 T *c,
323 const int &ldc);
324};
325
326template <typename T, typename Device>
328 using Real = typename GetTypeReal<T>::type;
329
330 void operator()(const Device *d, const int &nbase, const int &nbase_x, const int &notconv,
331 T *result, const T *vectors, const Real *eigenvalues);
332};
333
334template <typename T, typename Device>
336 using Real = typename GetTypeReal<T>::type;
337 void operator()(const Device* d,
338 const int& dim,
339 T* psi_iter,
340 const int& nbase,
341 const int& notconv,
342 const Real* precondition,
343 const Real* eigenvalues);
344};
345
346template <typename T, typename Device>
348 using Real = typename GetTypeReal<T>::type;
349 void operator()(const Device* d,
350 const int& dim,
351 T* psi_iter,
352 const int& nbase,
353 const int& notconv,
354 Real* psi_norm = nullptr);
355};
356
357template <typename T>
358struct normalize_op<T, base_device::DEVICE_GPU> {
359 using Real = typename GetTypeReal<T>::type;
360 void operator()(const base_device::DEVICE_GPU* d,
361 const int& dim,
362 T* psi_iter,
363 const int& nbase,
364 const int& notconv,
365 Real* psi_norm);
366};
367
368#if __CUDA || __UT_USE_CUDA || __ROCM || __UT_USE_ROCM
369// Partially specialize functor for base_device::GpuDevice.
370template <typename T> struct dot_real_op<T, base_device::DEVICE_GPU> {
371 using Real = typename GetTypeReal<T>::type;
372 Real operator()(const int &dim,
373 const T *psi_L, const T *psi_R, const bool reduce = true);
374};
375
376// vector operator: result[i] = vector[i] / constant
377template <typename T>
378struct vector_mul_real_op<T, base_device::DEVICE_GPU>
379{
380 using Real = typename GetTypeReal<T>::type;
381 void operator()(const int dim, T* result, const T* vector, const Real constant);
382};
383
384// vector operator: result[i] = vector1[i](complex) * vector2[i](not complex)
385template <typename T> struct vector_mul_vector_op<T, base_device::DEVICE_GPU> {
386 using Real = typename GetTypeReal<T>::type;
387 void operator()(const int& dim, T* result, const T* vector1, const Real* vector2, const bool& add = false);
388};
389
390// vector operator: result[i] = vector[i] / constant
391template <typename T> struct vector_div_constant_op<T, base_device::DEVICE_GPU> {
392 using Real = typename GetTypeReal<T>::type;
393 void operator()(const int& dim, T* result, const T* vector, const Real constant);
394};
395
396// vector operator: result[i] = vector1[i](complex) / vector2[i](not complex)
397template <typename T> struct vector_div_vector_op<T, base_device::DEVICE_GPU> {
398 using Real = typename GetTypeReal<T>::type;
399 void operator()(const int &dim, T *result,
400 const T *vector1, const Real *vector2);
401};
402
403// vector operator: result[i] = vector1[i] * constant1 + vector2[i] * constant2
404template <typename T>
405struct vector_add_vector_op<T, base_device::DEVICE_GPU> {
406 using Real = typename GetTypeReal<T>::type;
407 void operator()(const int &dim, T *result,
408 const T *vector1, const Real constant1, const T *vector2,
409 const Real constant2);
410};
411
412template <typename T> struct matrixCopy<T, base_device::DEVICE_GPU> {
413 void operator()(const int& n1,
414 const int& n2,
415 const T* A, // input
416 const int& LDA,
417 T* B, // output
418 const int& LDB);
419};
420
421template <typename T> struct matrix_mul_vector_op<T, base_device::DEVICE_GPU> {
422 using Real = typename GetTypeReal<T>::type;
423 void operator()(const int &m, const int &n,
424 T *a,
425 const int &lda,
426 const Real *b,
427 const Real alpha,
428 T *c,
429 const int &ldc);
430};
431
432void createGpuBlasHandle();
433void destoryBLAShandle();
434
435// vector operator: result[i] = -lambda[i] * vector[i]
436template <typename T> struct apply_eigenvalues_op<T, base_device::DEVICE_GPU> {
437 using Real = typename GetTypeReal<T>::type;
438
439 void operator()(const base_device::DEVICE_GPU *d, const int &nbase, const int &nbase_x, const int &notconv,
440 T *result, const T *vectors, const Real *eigenvalues);
441};
442
443template <typename T>
444struct precondition_op<T, base_device::DEVICE_GPU> {
445 using Real = typename GetTypeReal<T>::type;
446 void operator()(const base_device::DEVICE_GPU* d,
447 const int& dim,
448 T* psi_iter,
449 const int& nbase,
450 const int& notconv,
451 const Real* precondition,
452 const Real* eigenvalues);
453};
454
455#endif // __CUDA || __UT_USE_CUDA || __ROCM || __UT_USE_ROCM
456} // namespace hsolver
457
458#endif // MODULE_HSOLVER_MATH_KERNEL_H
#define N
Definition exp.cpp:24
#define T
Definition exp.cpp:237
Definition array_pool.h:6
std::complex< double > set_real_tocomplex(const std::complex< double > &x)
Definition math_kernel_op.h:24
std::complex< double > get_conj(const std::complex< double > &x)
Definition math_kernel_op.h:36
Definition device.cpp:21
T type
Definition macros.h:8
Definition math_kernel_op.h:327
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:328
Definition math_kernel_op.h:132
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:168
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:169
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:217
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:192
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:286
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:272
void operator()(const int &row, const int &col, const T *input_matrix, T *output_matrix)
transpose the input matrix
Definition math_kernel_op.h:302
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:303
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:359
Definition math_kernel_op.h:347
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:348
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:335
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:336
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:52
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:151
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:152
Definition math_kernel_op.h:101
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:102
void operator()(const int &dim, T *result, const T *vector, const Real constant)
result[i] = vector[i] / constant
Definition math_kernel_op.h:116
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:117
Definition math_kernel_op.h:68
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:69
Definition math_kernel_op.h:85
typename GetTypeReal< T >::type Real
Definition math_kernel_op.h:86
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)