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/realarray.h"
10
11#include "../../source_cell/unitcell.h"
12
13//==========================================================
14// CLASS :
15// NAME : Bessel_Basis
16//==========================================================
18{
19public:
22
37 void init(
38 const bool start_from_file,
39 const double &ecutwfc,
40 const int &ntype,
41 const int &lmax_in,
42 const bool &smooth,
43 const double &sigma,
44 const double &rcut_in,
45 const double &tol_in,
46 const UnitCell& ucell,
47 const double &dk = 0.01,
48 const double &dr = 0.01
49 );
53 const int& get_ecut_number() const { return Ecut_number;}
54
61 double Polynomial_Interpolation(const int &it, const int &l, const int &ic, const double &gnorm)const;
67 double Polynomial_Interpolation2(const int &l, const int &ie, const double &gnorm)const;
68
69
73 const double &get_ecut() const {return ecut;}
77 const double &get_rcut() const {return rcut;}
78
79 const double &get_tolerence() const {return tolerence;}
80
81
86 const bool &get_smooth() const {return smooth;}
90 const double &get_sigma() const {return sigma;}
91
92private:
95
98
101
103 int kmesh=0;
105 double Dk;
109 double rcut=0.0;
111 double ecut=0.0;
112 double tolerence=0.0;
114 bool smooth=false;
116 double sigma=0.0;
117
123 void allocate_C4(
124 const int &ntype,
125 const int &lmax,
126 const int &nmax,
127 const int &ecut_number,
128 const UnitCell& ucell
129 );
130
138 void readin_C4(
139 const std::string &name,
140 const int &ntype,
141 const int &ecut,
142 const int &rcut,
143 const int &ecut_number,
144 const double &tolerence,
145 const UnitCell& ucell
146 );
147
149
155 void init_Faln(
156 const int &ntype,
157 const int &lmax,
158 const int &nmax,
159 const int &ecut_number,
160 const UnitCell& ucell
161 );
162
164 int nwfc=0;
165
177 void init_TableOne(
178 const bool smooth_in,
179 const double &sigma_in,
180 const double &ecut,
181 const double &rcut,
182 const double &dr,
183 const double &dk,
184 const int &lmax,
185 const int &ecut_number,
186 const double &tolerence);
187};
188
189#endif
Definition bessel_basis.h:18
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:349
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:108
const double & get_sigma() const
get sigma the stddev (standard deviation) used in smooth function (Gaussian function)
Definition bessel_basis.h:90
const double & get_ecut() const
get energy cutoff, which is used to truncate SBF Jlq.
Definition bessel_basis.h:73
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:22
ModuleBase::realArray C4
Coefficients to be optimized!
Definition bessel_basis.h:97
ModuleBase::realArray Faln
the most important array to calculate spillage, has dimension (ntype, lmax+1, max_n,...
Definition bessel_basis.h:94
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:100
int kmesh
mesh of k vector, k is in j_l(k*r)
Definition bessel_basis.h:103
double Dk
grid of k
Definition bessel_basis.h:105
double Polynomial_Interpolation2(const int &l, const int &ie, const double &gnorm) const
Cubic spline interpolation for matrix TableOne.
Definition bessel_basis.cpp:83
const int & get_ecut_number() const
return number of SBFs used for one chi (see details for more information)
Definition bessel_basis.h:53
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:109
const double & get_rcut() const
cutoff radius of radial SBF Jlq.
Definition bessel_basis.h:77
Bessel_Basis()
Definition bessel_basis.cpp:10
~Bessel_Basis()
Definition bessel_basis.cpp:16
int Ecut_number
number of q vector, q is in j_l(q*r)
Definition bessel_basis.h:107
int nwfc
number of localized wave functions
Definition bessel_basis.h:164
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:461
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:126
double sigma
stddev of smooth function (Gaussian function, centered at rcut)
Definition bessel_basis.h:116
const bool & get_smooth() const
check if SBFs are smoothed (mohan add 2009-08-28)
Definition bessel_basis.h:86
double ecut
energy cutoff for determining kmesh and number of SBFs
Definition bessel_basis.h:111
bool smooth
whether smooth SBFs around cutoff radius, resulting in non-zero values. For importance of smooth of S...
Definition bessel_basis.h:114
double tolerence
Definition bessel_basis.h:112
const double & get_tolerence() const
Definition bessel_basis.h:79
double float array
Definition realarray.h:21
Definition unitcell.h:15