ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
radial_proj.h
Go to the documentation of this file.
1#ifndef RADIAL_PROJECTION_H
2#define RADIAL_PROJECTION_H
3
15#include "source_base/vector3.h"
17#include <memory>
18#include <vector>
19#include <complex>
20#include <map>
21#include <utility>
22#include <algorithm>
23
26
28{
66 {
67 public:
117 static void _build_backward_map(const std::vector<std::vector<int>>& it2iproj,
118 const std::vector<int>& iproj2l,
119 std::vector<int>& irow2it,
120 std::vector<int>& irow2iproj,
121 std::vector<int>& irow2m);
122 static void _build_forward_map(const std::vector<std::vector<int>>& it2ia,
123 const std::vector<std::vector<int>>& it2iproj,
124 const std::vector<int>& iproj2l,
125 std::map<std::tuple<int, int, int, int>, int>& itiaiprojm2irow);
126 public:
135 RadialProjector(const bool realspace = false) {}
137
138 // it is more feasible to build interpolation table. this function will tabulate
139 // for those functions. This function will write the in-build tab_, to place the
140 // values of Jl[f](q) for each q and l.
141 // Here I provide two versions of tabulate, one for the case may be capable to
142 // avoid the memory copy operation, and the other one is for the case of
143 // std::vector<double> and std::vector<std::vector<double>>.
154 void _build_sbt_tab(const int nr,
155 const double* r,
156 const std::vector<double*>& radials,
157 const std::vector<int>& l,
158 const int nq, //< GlobalV::DQ
159 const double& dq); //< GlobalV::NQX
160 void _build_sbt_tab(const std::vector<double>& r,
161 const std::vector<std::vector<double>>& radials,
162 const std::vector<int>& l,
163 const int nq, //< GlobalV::DQ
164 const double& dq); //< GlobalV::NQX
165 // compatibility concern: for FS_Nonlocal_tools. Will not call sbtft so need omega
166 void _build_sbt_tab(const std::vector<int>& nproj,
167 const std::vector<double>& r,
168 const std::vector<std::vector<double>>& radials,
169 const std::vector<int>& l,
170 const int nq, //< GlobalV::DQ
171 const double& dq, //< GlobalV::NQX
172 const double& omega,
173 const int npol,
175 ModuleBase::matrix& nhtol);
188 void sbtft(const std::vector<ModuleBase::Vector3<double>>& qs,
189 std::vector<std::complex<double>>& out,
190 const char type = 'r', // 'r' for ket |>, 'l' for bra <|
191 const double& omega = 1.0,
192 const double& tpiba = 1.0); // 'n' for no gradient, 'x', 'y', 'z' for gradient in x, y, z direction
193
194 void sbfft(); // interface for SBFFT
195
196
197 private:
198 std::unique_ptr<ModuleBase::CubicSpline> cubspl_;
199 std::vector<int> l_;
200 };
201
232 void _mask_func(std::vector<double>& mask);
233
246 void _do_mask_on_radial(const int nr1,
247 const double* r,
248 const double* in,
249 const int nr2,
250 const double* mask,
251 double* out);
252}
253
254#endif // RADIAL_PROJECTION_H
3 elements vector
Definition vector3.h:24
Definition matrix.h:18
double float array
Definition realarray.h:21
RadialProjector is for projecting a function who has seperatable radial and angular parts: f(r) = f(|...
Definition radial_proj.h:66
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)
Definition radial_proj.cpp:53
~RadialProjector()
Definition radial_proj.h:136
RadialProjector(const bool realspace=false)
Construct a new Radial Projector object.
Definition radial_proj.h:135
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)
Definition radial_proj.cpp:78
std::unique_ptr< ModuleBase::CubicSpline > cubspl_
Definition radial_proj.h:198
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(ome...
Definition radial_proj.cpp:195
std::vector< int > l_
Definition radial_proj.h:199
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)
Definition radial_proj.cpp:13
Definition radial_proj.h:28
void _mask_func(std::vector< double > &mask)
get the mask function for SBFFT
Definition radial_proj.cpp:243
void _do_mask_on_radial(const int nr1, const double *r, const double *in, const int nr2, const double *mask, double *out)
do operation w(r)/m(r) on a radial function. The cutoff radius of w(r) is smaller than the cutoff rad...
Definition radial_proj.cpp:305