ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
|
RadialProjector is for projecting a function who has seperatable radial and angular parts: f(r) = f(|r|) * Ylm(theta, phi) onto planewave basis. More...
#include <radial_proj.h>
Public Member Functions | |
RadialProjector (const bool realspace=false) | |
Construct a new Radial Projector object. | |
~RadialProjector () | |
void | _build_sbt_tab (const int nr, const double *r, const std::vector< double * > &radials, const std::vector< int > &l, const int nq, const double &dq) |
make a interpolation table for the Spherical Bessel Transform of f(r) | |
void | _build_sbt_tab (const std::vector< double > &r, const std::vector< std::vector< double > > &radials, const std::vector< int > &l, const int nq, const double &dq) |
void | _build_sbt_tab (const std::vector< int > &nproj, const std::vector< double > &r, const std::vector< std::vector< double > > &radials, const std::vector< int > &l, const int nq, const double &dq, const double &omega, const int npol, ModuleBase::realArray &tab, ModuleBase::matrix &nhtol) |
void | sbtft (const std::vector< ModuleBase::Vector3< double > > &qs, std::vector< std::complex< double > > &out, const char type='r', const double &omega=1.0, const double &tpiba=1.0) |
perform analytical version of the Fourier transform: F(q) = int(f(r)*exp(-iq.r) d^3r) = 4*pi/sqrt(omega) * (-i)^l * Jl[f](q) * Ylm(q) , where Ylm(q) is real spherical harmonic function, and Jl[f](q) is the Spherial Bessel Transform of f(r): Jl[f](q) = int(f(r)*j_l(q*r)*r^2 dr) , where j_l(q*r) is the spherical Bessel function of the first kind. . If use another notation, F(q) = <q|f>, this is denoted as type "r" for ket |>, and "l" for bra <|. | |
void | sbfft () |
Static Public Member Functions | |
static void | _build_backward_map (const std::vector< std::vector< int > > &it2iproj, const std::vector< int > &iproj2l, std::vector< int > &irow2it, std::vector< int > &irow2iproj, std::vector< int > &irow2m) |
static void | _build_forward_map (const std::vector< std::vector< int > > &it2ia, const std::vector< std::vector< int > > &it2iproj, const std::vector< int > &iproj2l, std::map< std::tuple< int, int, int, int >, int > &itiaiprojm2irow) |
Private Attributes | |
std::unique_ptr< ModuleBase::CubicSpline > | cubspl_ |
std::vector< int > | l_ |
RadialProjector is for projecting a function who has seperatable radial and angular parts: f(r) = f(|r|) * Ylm(theta, phi) onto planewave basis.
Usage:
|
inline |
Construct a new Radial Projector object.
realspace | if perform integration in real space rather than reciprocal space , default is false |
|
inline |
|
static |
Notation of following two functions:
Given all the projectors are listed in a series, so the iproj
is the index goes across all atomtypes, which means if for the first type, the iproj goes from 0 to 4, then the second atomtypes the iproj will start from 5, and so on... However, there is also another convention, like numerical atomic orbitals, developer always use "l" to index orbitals, here, in all output map, the iproj
will start from 0, which
First, the following lists should be prepared as early as possible,
it2iproj: for given it, the index of atom type, return the list of index of projectors.
iproj2l: for given iproj, the index of projectors, return the l of this projector. More simply explaning, it is just the list of angular momentum of projectors.
it2ia: just a list that stolen information from UnitCell, for given it, the index of atom within the range of it. So this list is different from the it2iproj, iproj is the index across type but ia is the index within the type. So for each it2ia[it], the ia, in principle , always/can start from 0.
One may question that does the indexing support one atom type with multiple projectors? The answer is YES. Combining the it2iproj and it2ia, one can even support PART of atoms of one
Then the returned lists,
irow2it: for given irow
, the index of row, return the it
: the index of atom type.
irow2ia: for given irow
, the index of row, return the ia
: the index of atom within the range of it
.
irow2iproj: for given irow
, the index of row, return the iproj
, the index of projectors, note that this iproj
is the local index.
irow2m: for given irow, the index of row, return the m, the magnetic quantum number of this projector.
One may complain that cannot get l
from the irow
, but the truth is, not exactly. One can get the l
starting from irow
by:
|
static |
void RadialProjection::RadialProjector::_build_sbt_tab | ( | const int | nr, |
const double * | r, | ||
const std::vector< double * > & | radials, | ||
const std::vector< int > & | l, | ||
const int | nq, | ||
const double & | dq | ||
) |
make a interpolation table for the Spherical Bessel Transform of f(r)
nr | number of grid points, shared by all radial functions |
r | radial grids, shared by all radial functions |
radials | radial functions, each element is a radial function |
l | angular momentum quantum number for each radial function |
nq | number of q-points |
dq | space between q-points |
void RadialProjection::RadialProjector::_build_sbt_tab | ( | const std::vector< double > & | r, |
const std::vector< std::vector< double > > & | radials, | ||
const std::vector< int > & | l, | ||
const int | nq, | ||
const double & | dq | ||
) |
void RadialProjection::RadialProjector::_build_sbt_tab | ( | const std::vector< int > & | nproj, |
const std::vector< double > & | r, | ||
const std::vector< std::vector< double > > & | radials, | ||
const std::vector< int > & | l, | ||
const int | nq, | ||
const double & | dq, | ||
const double & | omega, | ||
const int | npol, | ||
ModuleBase::realArray & | tab, | ||
ModuleBase::matrix & | nhtol | ||
) |
void RadialProjection::RadialProjector::sbfft | ( | ) |
void RadialProjection::RadialProjector::sbtft | ( | const std::vector< ModuleBase::Vector3< double > > & | qs, |
std::vector< std::complex< double > > & | out, | ||
const char | type = 'r' , |
||
const double & | omega = 1.0 , |
||
const double & | tpiba = 1.0 |
||
) |
perform analytical version of the Fourier transform: F(q) = int(f(r)*exp(-iq.r) d^3r) = 4*pi/sqrt(omega) * (-i)^l * Jl[f](q) * Ylm(q) , where Ylm(q) is real spherical harmonic function, and Jl[f](q) is the Spherial Bessel Transform of f(r): Jl[f](q) = int(f(r)*j_l(q*r)*r^2 dr) , where j_l(q*r) is the spherical Bessel function of the first kind. . If use another notation, F(q) = <q|f>, this is denoted as type "r" for ket |>, and "l" for bra <|.
|
private |
|
private |