ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
pw_basis_k.h
Go to the documentation of this file.
1#ifndef PWBASISK_H
2#define PWBASISK_H
3
4#include "pw_basis.h"
6namespace ModulePW
7{
8
55class PW_Basis_K : public PW_Basis
56{
57
58public:
59 PW_Basis_K();
60 PW_Basis_K(std::string device_, std::string precision_) : PW_Basis(device_, precision_) {classname="PW_Basis_K";}
62
63 //init parameters of pw_basis_k class
64 void initparameters(
65 const bool gamma_only_in,
66 const double ecut_in,
67 const int nk_in, //number of k points in this pool
68 const ModuleBase::Vector3<double> *kvec_d, // Direct coordinates of k points
69 const int distribution_type_in = 1,
70 const bool xprime_in = true
71 );
72
73 public:
74 int nks=0;//number of k points in this pool
75 ModuleBase::Vector3<double> *kvec_d=nullptr; // Direct coordinates of k points
76 ModuleBase::Vector3<double> *kvec_c=nullptr; // Cartesian coordinates of k points
77 int *npwk=nullptr; //[nks] number of plane waves of different k-points
78 int npwk_max=0; //max npwk among all nks k-points, it may be smaller than npw
79 //npw cutoff: (|g|+|k|)^2, npwk in the the npw ball, thus is smaller
80 double gk_ecut=0; //Energy cut off for (g+k)^2/2
81
82public:
83 //prepare for transforms between real and reciprocal spaces
84 void setuptransform();
85
86 int *igl2isz_k=nullptr, * d_igl2isz_k = nullptr; //[npwk_max*nks] map (igl,ik) to (is,iz)
87 int *igl2ig_k=nullptr;//[npwk_max*nks] map (igl,ik) to ig
88 int *ig2ixyz_k=nullptr;
89 std::vector<int> ig2ixyz_k_cpu;
90 double *gk2=nullptr; // modulus (G+K)^2 of G vectors [npwk_max*nks]
91
92 // liuyu add 2023-09-06
93 double erf_ecut=0.0; // the value of the constant energy cutoff
94 double erf_height=0.0; // the height of the energy step for reciprocal vectors
95 double erf_sigma=0.0; // the width of the energy step for reciprocal vectors
96
97 //collect gdirect, gcar, gg
98 void collect_local_pw(const double& erf_ecut_in = 0.0,
99 const double& erf_height_in = 0.0,
100 const double& erf_sigma_in = 0.1);
101
102 private:
103 float * s_gk2 = nullptr;
104 double * d_gk2 = nullptr; // modulus (G+K)^2 of G vectors [npwk_max*nks]
105 //create igl2isz_k map array for fft
106 void setupIndGk();
107 // get ig2ixyz_k
108 void get_ig2ixyz_k();
109 //calculate G+K, it is a private function
110 ModuleBase::Vector3<double> cal_GplusK_cartesian(const int ik, const int ig) const;
111
112 public:
113 template <typename FPTYPE>
114 void real2recip(const FPTYPE* in,
115 std::complex<FPTYPE>* out,
116 const int ik,
117 const bool add = false,
118 const FPTYPE factor = 1.0) const; // in:(nplane,nx*ny) ; out(nz, ns)
119 template <typename FPTYPE>
120 void real2recip(const std::complex<FPTYPE>* in,
121 std::complex<FPTYPE>* out,
122 const int ik,
123 const bool add = false,
124 const FPTYPE factor = 1.0) const; // in:(nplane,nx*ny) ; out(nz, ns)
125 template <typename FPTYPE>
126 void recip2real(const std::complex<FPTYPE>* in,
127 FPTYPE* out,
128 const int ik,
129 const bool add = false,
130 const FPTYPE factor = 1.0) const; // in:(nz, ns) ; out(nplane,nx*ny)
131 template <typename FPTYPE>
132 void recip2real(const std::complex<FPTYPE>* in,
133 std::complex<FPTYPE>* out,
134 const int ik,
135 const bool add = false,
136 const FPTYPE factor = 1.0) const; // in:(nz, ns) ; out(nplane,nx*ny)
137 #if defined(__DSP)
138 template <typename FPTYPE, typename Device>
139 void convolution(const Device* ctx,
140 const int ik,
141 const int size,
142 const std::complex<FPTYPE>* input,
143 const FPTYPE* input1,
144 std::complex<FPTYPE>* output,
145 const bool add = false,
146 const FPTYPE factor =1.0) const ;
147
148 template <typename FPTYPE>
149 void real2recip_dsp(const std::complex<FPTYPE>* in,
150 std::complex<FPTYPE>* out,
151 const int ik,
152 const bool add = false,
153 const FPTYPE factor = 1.0) const; // in:(nplane,nx*ny) ; out(nz, ns)
154 template <typename FPTYPE>
155 void recip2real_dsp(const std::complex<FPTYPE>* in,
156 std::complex<FPTYPE>* out,
157 const int ik,
158 const bool add = false,
159 const FPTYPE factor = 1.0) const; // in:(nz, ns) ; out(nplane,nx*ny)
160
161 #endif
162
163 template <typename FPTYPE, typename Device>
164 void real_to_recip(const Device* ctx,
165 const std::complex<FPTYPE>* in,
166 std::complex<FPTYPE>* out,
167 const int ik,
168 const bool add = false,
169 const FPTYPE factor = 1.0) const; // in:(nplane,nx*ny) ; out(nz, ns)
170 template <typename FPTYPE, typename Device>
171 void recip_to_real(const Device* ctx,
172 const std::complex<FPTYPE>* in,
173 std::complex<FPTYPE>* out,
174 const int ik,
175 const bool add = false,
176 const FPTYPE factor = 1.0) const; // in:(nz, ns) ; out(nplane,nx*ny)
177
178
179 template <typename TK,
180 typename Device,
181 typename std::enable_if<std::is_same<Device, base_device::DEVICE_CPU>::value, int>::type = 0>
182 void real_to_recip(const TK* in,
183 TK* out,
184 const int ik,
185 const bool add = false,
186 const typename GetTypeReal<TK>::type factor = 1.0) const
187 {
188 #if defined(__DSP)
189 this->real2recip_dsp(in, out, ik, add, factor);
190 #else
191 this->real2recip(in,out,ik,add,factor);
192 #endif
193 }
194 template <typename TK,
195 typename Device,
196 typename std::enable_if<std::is_same<Device, base_device::DEVICE_CPU>::value, int>::type = 0>
197 void recip_to_real(const TK* in,
198 TK* out,
199 const int ik,
200 const bool add = false,
201 const typename GetTypeReal<TK>::type factor = 1.0) const
202 {
203
204 #if defined(__DSP)
205 this->recip2real_dsp(in,out,ik,add,factor);
206 #else
207 this->recip2real(in,out,ik,add,factor);
208 #endif
209 }
210 template <typename FPTYPE>
211 void real2recip_gpu(const std::complex<FPTYPE>* in,
212 std::complex<FPTYPE>* out,
213 const int ik,
214 const bool add = false,
215 const FPTYPE factor = 1.0) const; // in:(nplane,nx*ny) ; out(nz, ns)
216
217 template <typename FPTYPE>
218 void recip2real_gpu(const std::complex<FPTYPE>* in,
219 std::complex<FPTYPE>* out,
220 const int ik,
221 const bool add = false,
222 const FPTYPE factor = 1.0) const; // in:(nz, ns) ; out(nplane,nx*ny)
223
224 template <typename FPTYPE,
225 typename Device,
226 typename std::enable_if<!std::is_same<Device, base_device::DEVICE_CPU>::value, int>::type = 0>
227 void real_to_recip(const FPTYPE* in,
228 FPTYPE* out,
229 const int ik,
230 const bool add = false,
231 const typename GetTypeReal<FPTYPE>::type factor = 1.0) const
232 {
233 this->real2recip_gpu(in, out, ik, add, factor);
234 }
235
236 template <typename TK,
237 typename Device,
238 typename std::enable_if<std::is_same<Device, base_device::DEVICE_GPU>::value, int>::type = 0>
239 void recip_to_real(const TK* in,
240 TK* out,
241 const int ik,
242 const bool add = false,
243 const typename GetTypeReal<TK>::type factor = 1.0) const
244 {
245 this->recip2real_gpu(in, out, ik, add, factor);
246 }
247
248 public:
249 //operator:
250 //get (G+K)^2:
251 double& getgk2(const int ik, const int igl) const;
252 //get G
253 ModuleBase::Vector3<double>& getgcar(const int ik, const int igl) const;
254 //get G-direct
255 ModuleBase::Vector3<double> getgdirect(const int ik, const int igl) const;
256 //get (G+K)
257 ModuleBase::Vector3<double> getgpluskcar(const int ik, const int igl) const;
258 //get igl2isz_k
259 int& getigl2isz(const int ik, const int igl) const;
260 //get igl2ig_k or igk(ik,ig) in older ABACUS
261 int& getigl2ig(const int ik, const int igl) const;
262
263 //get ig_to_ix
264 std::vector<int> get_ig2ix(const int ik) const;
265 //get ig_to_iy
266 std::vector<int> get_ig2iy(const int ik) const;
267 //get ig_to_iz
268 std::vector<int> get_ig2iz(const int ik) const;
269
270 template <typename FPTYPE> FPTYPE * get_gk2_data() const;
271 template <typename FPTYPE> FPTYPE * get_gcar_data() const;
272 template <typename FPTYPE> FPTYPE * get_kvec_c_data() const;
273
274private:
275 float * s_gcar = nullptr, * s_kvec_c = nullptr;
276 double * d_gcar = nullptr, * d_kvec_c = nullptr;
277};
278
279}
280#endif //PlaneWave_K class
281
282#include "./pw_basis_k_big.h" //temporary it will be removed
283
3 elements vector
Definition vector3.h:24
Special pw_basis class. It includes different k-points.
Definition pw_basis_k.h:56
double * d_kvec_c
Definition pw_basis_k.h:276
ModuleBase::Vector3< double > getgdirect(const int ik, const int igl) const
Definition pw_basis_k.cpp:390
float * s_gk2
Definition pw_basis_k.h:103
FPTYPE * get_kvec_c_data() const
void real_to_recip(const Device *ctx, const std::complex< FPTYPE > *in, std::complex< FPTYPE > *out, const int ik, const bool add=false, const FPTYPE factor=1.0) const
Definition xc3_mock.h:81
ModuleBase::Vector3< double > & getgcar(const int ik, const int igl) const
Definition pw_basis_k.cpp:385
PW_Basis_K()
Definition pw_basis_k.cpp:12
int npwk_max
Definition pw_basis_k.h:78
int * ig2ixyz_k
[npw] map ig to ixyz
Definition pw_basis_k.h:88
void recip2real(const std::complex< FPTYPE > *in, FPTYPE *out, const int ik, const bool add=false, const FPTYPE factor=1.0) const
transform reciprocal space to real space
Definition pw_transform_k.cpp:228
~PW_Basis_K()
Definition pw_basis_k.cpp:17
void real_to_recip(const FPTYPE *in, FPTYPE *out, const int ik, const bool add=false, const typename GetTypeReal< FPTYPE >::type factor=1.0) const
Definition pw_basis_k.h:227
double gk_ecut
Definition pw_basis_k.h:80
std::vector< int > get_ig2ix(const int ik) const
Definition pw_basis_k.cpp:441
double * d_gk2
Definition pw_basis_k.h:104
double & getgk2(const int ik, const int igl) const
Definition pw_basis_k.cpp:380
void recip_to_real(const Device *ctx, const std::complex< FPTYPE > *in, std::complex< FPTYPE > *out, const int ik, const bool add=false, const FPTYPE factor=1.0) const
Definition xc3_mock.h:94
double erf_ecut
Definition pw_basis_k.h:93
void real2recip(const FPTYPE *in, std::complex< FPTYPE > *out, const int ik, const bool add=false, const FPTYPE factor=1.0) const
transform real space to reciprocal space
Definition pw_transform_k.cpp:91
int * npwk
Definition pw_basis_k.h:77
void setuptransform()
Definition pw_basis_k.cpp:206
std::vector< int > get_ig2iz(const int ik) const
Definition pw_basis_k.cpp:481
ModuleBase::Vector3< double > cal_GplusK_cartesian(const int ik, const int ig) const
Definition pw_basis_k.cpp:352
double * d_gcar
Definition pw_basis_k.h:276
void real2recip_gpu(const std::complex< FPTYPE > *in, std::complex< FPTYPE > *out, const int ik, const bool add=false, const FPTYPE factor=1.0) const
int & getigl2ig(const int ik, const int igl) const
Definition pw_basis_k.cpp:407
double erf_height
Definition pw_basis_k.h:94
int * igl2ig_k
Definition pw_basis_k.h:87
void initparameters(const bool gamma_only_in, const double ecut_in, const int nk_in, const ModuleBase::Vector3< double > *kvec_d, const int distribution_type_in=1, const bool xprime_in=true)
Definition pw_basis_k.cpp:49
double erf_sigma
Definition pw_basis_k.h:95
int * igl2isz_k
Definition pw_basis_k.h:86
void recip_to_real(const TK *in, TK *out, const int ik, const bool add=false, const typename GetTypeReal< TK >::type factor=1.0) const
Definition pw_basis_k.h:197
float * s_gcar
Definition pw_basis_k.h:275
double * gk2
[npw] map ig to ixyz,which is used in dsp fft.
Definition pw_basis_k.h:90
ModuleBase::Vector3< double > getgpluskcar(const int ik, const int igl) const
Definition pw_basis_k.cpp:399
std::vector< int > get_ig2iy(const int ik) const
Definition pw_basis_k.cpp:461
ModuleBase::Vector3< double > * kvec_d
Definition pw_basis_k.h:75
PW_Basis_K(std::string device_, std::string precision_)
Definition pw_basis_k.h:60
int * d_igl2isz_k
Definition pw_basis_k.h:86
std::vector< int > ig2ixyz_k_cpu
Definition pw_basis_k.h:89
int & getigl2isz(const int ik, const int igl) const
Definition pw_basis_k.cpp:403
ModuleBase::Vector3< double > * kvec_c
Definition pw_basis_k.h:76
float * s_kvec_c
Definition pw_basis_k.h:275
void get_ig2ixyz_k()
Definition pw_basis_k.cpp:412
FPTYPE * get_gk2_data() const
void recip2real_gpu(const std::complex< FPTYPE > *in, std::complex< FPTYPE > *out, const int ik, const bool add=false, const FPTYPE factor=1.0) const
void real_to_recip(const TK *in, TK *out, const int ik, const bool add=false, const typename GetTypeReal< TK >::type factor=1.0) const
Definition pw_basis_k.h:182
void setupIndGk()
Definition pw_basis_k.cpp:130
FPTYPE * get_gcar_data() const
int nks
Definition pw_basis_k.h:74
A class which can convert a function of "r" to the corresponding linear superposition of plane waves ...
Definition pw_basis.h:56
std::string classname
Definition pw_basis.h:59
void collect_local_pw()
Definition pw_basis.cpp:135
Definition output.h:13
Definition pw_op.cpp:3
T type
Definition macros.h:8