ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
structure_factor.h
Go to the documentation of this file.
1#ifndef STRUCTURE_FACTOR_H
2#define STRUCTURE_FACTOR_H
3
8#include "source_psi/psi.h"
9
11{
12
13public:
16 void set(const ModulePW::PW_Basis* rho_basis_in, const int& nbspline_in);
17
18 //===============================================
19 // Part 4: G vectors in reciprocal FFT box
20 //===============================================
21 public:
22 int nbspline=0;
23
24 // structure factor (ntype, ngmc)
26
27 void setup(const UnitCell* Ucell,
28 const Parallel_Grid& pgrid,
29 const ModulePW::PW_Basis* rho_basis); // Calculate structure factors
30
32 void bspline_sf(
33 const int,
34 const UnitCell* Ucell,
35 const Parallel_Grid& pgrid,
37
38 void bsplinecoef(std::complex<double> *b1, std::complex<double> *b2, std::complex<double> *b3,
39 const int nx, const int ny, const int nz, const int norder);
40
41
42public:
43 // phase of e^{-iG*tau_s}
44 ModuleBase::ComplexMatrix eigts1; // dimension: [Ucell->nat, 2*this->ncx + 1]
45 ModuleBase::ComplexMatrix eigts2; // dimension: [Ucell->nat, 2*this->ncy + 1]
46 ModuleBase::ComplexMatrix eigts3; // dimension: [Ucell->nat, 2*this->ncz + 1]
47
48 template <typename FPTYPE> std::complex<FPTYPE> * get_eigts1_data() const;
49 template <typename FPTYPE> std::complex<FPTYPE> * get_eigts2_data() const;
50 template <typename FPTYPE> std::complex<FPTYPE> * get_eigts3_data() const;
51
52 public:
53 // sf with k points
54 std::complex<double>* get_sk(const int ik, const int it, const int ia, const ModulePW::PW_Basis_K* wfc_basis) const;
55 template <typename FPTYPE, typename Device>
56
57 void get_sk(Device* ctx, const int ik, const ModulePW::PW_Basis_K* wfc_basis, std::complex<FPTYPE>* sk) const;
58
59 std::complex<double>* get_skq(int ik,
60 int it,
61 int ia,
62 const ModulePW::PW_Basis_K* wfc_basis,
64
65 private:
66
67 const UnitCell* ucell=nullptr;
68 std::complex<float> * c_eigts1 = nullptr;
69 std::complex<float> * c_eigts2 = nullptr;
70 std::complex<float> * c_eigts3 = nullptr;
71
72 std::complex<double> * z_eigts1 = nullptr;
73 std::complex<double> * z_eigts2 = nullptr;
74 std::complex<double> * z_eigts3 = nullptr;
75
76 const ModulePW::PW_Basis* rho_basis = nullptr;
77 std::string device = "cpu";
78};
79#endif //PlaneWave class
Definition complexmatrix.h:14
3 elements vector
Definition vector3.h:22
Special pw_basis class. It includes different k-points.
Definition pw_basis_k.h:57
A class which can convert a function of "r" to the corresponding linear superposition of plane waves ...
Definition pw_basis.h:56
Definition parallel_grid.h:8
Definition structure_factor.h:11
ModuleBase::ComplexMatrix eigts3
Definition structure_factor.h:46
std::complex< double > * z_eigts3
Definition structure_factor.h:74
const UnitCell * ucell
Definition structure_factor.h:67
~Structure_Factor()
Definition charge_extra_test.cpp:89
std::complex< double > * z_eigts1
Definition structure_factor.h:72
void bsplinecoef(std::complex< double > *b1, std::complex< double > *b2, std::complex< double > *b3, const int nx, const int ny, const int nz, const int norder)
Definition structure_factor.cpp:311
std::complex< FPTYPE > * get_eigts3_data() const
void setup(const UnitCell *Ucell, const Parallel_Grid &pgrid, const ModulePW::PW_Basis *rho_basis)
Definition charge_extra_test.cpp:92
std::string device
Definition structure_factor.h:77
Structure_Factor()
Definition charge_extra_test.cpp:86
std::complex< float > * c_eigts2
Definition structure_factor.h:69
const ModulePW::PW_Basis * rho_basis
Definition structure_factor.h:76
std::complex< float > * c_eigts3
Definition structure_factor.h:70
std::complex< double > * get_skq(int ik, int it, int ia, const ModulePW::PW_Basis_K *wfc_basis, ModuleBase::Vector3< double > q)
Definition structure_factor_k.cpp:146
void set(const ModulePW::PW_Basis *rho_basis_in, const int &nbspline_in)
Definition structure_factor.cpp:47
std::complex< FPTYPE > * get_eigts1_data() const
std::complex< FPTYPE > * get_eigts2_data() const
int nbspline
Definition structure_factor.h:22
ModuleBase::ComplexMatrix eigts1
Definition structure_factor.h:44
std::complex< double > * get_sk(const int ik, const int it, const int ia, const ModulePW::PW_Basis_K *wfc_basis) const
Definition psi_initializer_unit_test.cpp:84
void bspline_sf(const int, const UnitCell *Ucell, const Parallel_Grid &pgrid, const ModulePW::PW_Basis *rho_basis)
calculate structure factors through Cardinal B-spline interpolation
Definition structure_factor.cpp:205
std::complex< float > * c_eigts1
Definition structure_factor.h:68
ModuleBase::ComplexMatrix eigts2
Definition structure_factor.h:45
std::complex< double > * z_eigts2
Definition structure_factor.h:73
ModuleBase::ComplexMatrix strucFac
Definition structure_factor.h:25
Definition unitcell.h:17