ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
bessel_basis.h
Go to the documentation of this file.
1//==========================================================
2// AUTHOR : mohan
3// DATE : 2009-3-29
4// Last Modify : 2021-01-04
5//==========================================================
6#ifndef BESSEL_BASIS_H
7#define BESSEL_BASIS_H
8#include "../source_base/global_function.h"
9#include "../source_base/global_variable.h"
10#include "../source_base/realarray.h"
11
12#include "../source_cell/unitcell.h"
13
14//==========================================================
15// CLASS :
16// NAME : Bessel_Basis
17//==========================================================
19{
20public:
23
38 void init(
39 const bool start_from_file,
40 const double &ecutwfc,
41 const int &ntype,
42 const int &lmax_in,
43 const bool &smooth,
44 const double &sigma,
45 const double &rcut_in,
46 const double &tol_in,
47 const UnitCell& ucell,
48 const double &dk = 0.01,
49 const double &dr = 0.01
50 );
54 const int& get_ecut_number() const { return Ecut_number;}
55
62 double Polynomial_Interpolation(const int &it, const int &l, const int &ic, const double &gnorm)const;
68 double Polynomial_Interpolation2(const int &l, const int &ie, const double &gnorm)const;
69
70
74 const double &get_ecut() const {return ecut;}
78 const double &get_rcut() const {return rcut;}
79
80 const double &get_tolerence() const {return tolerence;}
81
82
87 const bool &get_smooth() const {return smooth;}
91 const double &get_sigma() const {return sigma;}
92
93private:
96
99
102
104 int kmesh=0;
106 double Dk;
110 double rcut=0.0;
112 double ecut=0.0;
113 double tolerence=0.0;
115 bool smooth=false;
117 double sigma=0.0;
118
124 void allocate_C4(
125 const int &ntype,
126 const int &lmax,
127 const int &nmax,
128 const int &ecut_number,
129 const UnitCell& ucell
130 );
131
139 void readin_C4(
140 const std::string &name,
141 const int &ntype,
142 const int &ecut,
143 const int &rcut,
144 const int &ecut_number,
145 const double &tolerence,
146 const UnitCell& ucell
147 );
148
150
156 void init_Faln(
157 const int &ntype,
158 const int &lmax,
159 const int &nmax,
160 const int &ecut_number,
161 const UnitCell& ucell
162 );
163
165 int nwfc=0;
166
178 void init_TableOne(
179 const bool smooth_in,
180 const double &sigma_in,
181 const double &ecut,
182 const double &rcut,
183 const double &dr,
184 const double &dk,
185 const int &lmax,
186 const int &ecut_number,
187 const double &tolerence);
188};
189
190#endif
Definition bessel_basis.h:19
void readin_C4(const std::string &name, const int &ntype, const int &ecut, const int &rcut, const int &ecut_number, const double &tolerence, const UnitCell &ucell)
Read C4 from external file. Presently an O(N^2) search algorithm is used. A HTML parser is needed in ...
Definition bessel_basis.cpp:350
double Polynomial_Interpolation(const int &it, const int &l, const int &ic, const double &gnorm) const
Cubic spline interpolation for matrix Faln.
Definition bessel_basis.cpp:109
const double & get_sigma() const
get sigma the stddev (standard deviation) used in smooth function (Gaussian function)
Definition bessel_basis.h:91
const double & get_ecut() const
get energy cutoff, which is used to truncate SBF Jlq.
Definition bessel_basis.h:74
void init(const bool start_from_file, const double &ecutwfc, const int &ntype, const int &lmax_in, const bool &smooth, const double &sigma, const double &rcut_in, const double &tol_in, const UnitCell &ucell, const double &dk=0.01, const double &dr=0.01)
Initialization of Bessel function related matrices.
Definition bessel_basis.cpp:23
ModuleBase::realArray C4
Coefficients to be optimized!
Definition bessel_basis.h:98
ModuleBase::realArray Faln
the most important array to calculate spillage, has dimension (ntype, lmax+1, max_n,...
Definition bessel_basis.h:95
ModuleBase::realArray TableOne
matrix whose elements are int{dr r^2 j_l(qr)*j_l(kr)}, has dimension (lmax+1, nq, nk)
Definition bessel_basis.h:101
int kmesh
mesh of k vector, k is in j_l(k*r)
Definition bessel_basis.h:104
double Dk
grid of k
Definition bessel_basis.h:106
double Polynomial_Interpolation2(const int &l, const int &ie, const double &gnorm) const
Cubic spline interpolation for matrix TableOne.
Definition bessel_basis.cpp:84
const int & get_ecut_number() const
return number of SBFs used for one chi (see details for more information)
Definition bessel_basis.h:54
void init_TableOne()
double rcut
Cutoff radius (in a.u.) of SBFs, for any SBF j_l(qr), r>=rcut, j_l(q*r) = 0 (if not smoothed)
Definition bessel_basis.h:110
const double & get_rcut() const
cutoff radius of radial SBF Jlq.
Definition bessel_basis.h:78
Bessel_Basis()
Definition bessel_basis.cpp:11
~Bessel_Basis()
Definition bessel_basis.cpp:17
int Ecut_number
number of q vector, q is in j_l(q*r)
Definition bessel_basis.h:108
int nwfc
number of localized wave functions
Definition bessel_basis.h:165
void allocate_C4(const int &ntype, const int &lmax, const int &nmax, const int &ecut_number, const UnitCell &ucell)
Allocate memory for C4 matrix and initialize all elements to one.
Definition bessel_basis.cpp:462
void init_Faln(const int &ntype, const int &lmax, const int &nmax, const int &ecut_number, const UnitCell &ucell)
calculate F_{aln}(it, il, in, ik) = sum_{ie}{C4(it, il, in, ie)*TableOne(il, ie, ik)},...
Definition bessel_basis.cpp:127
double sigma
stddev of smooth function (Gaussian function, centered at rcut)
Definition bessel_basis.h:117
const bool & get_smooth() const
check if SBFs are smoothed (mohan add 2009-08-28)
Definition bessel_basis.h:87
double ecut
energy cutoff for determining kmesh and number of SBFs
Definition bessel_basis.h:112
bool smooth
whether smooth SBFs around cutoff radius, resulting in non-zero values. For importance of smooth of S...
Definition bessel_basis.h:115
double tolerence
Definition bessel_basis.h:113
const double & get_tolerence() const
Definition bessel_basis.h:80
double float array
Definition realarray.h:21
Definition unitcell.h:16