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
25#include "source_psi/psi.h"
27
29{
67 {
68 public:
118 static void _build_backward_map(const std::vector<std::vector<int>>& it2iproj,
119 const std::vector<int>& iproj2l,
120 std::vector<int>& irow2it,
121 std::vector<int>& irow2iproj,
122 std::vector<int>& irow2m);
123 static void _build_forward_map(const std::vector<std::vector<int>>& it2ia,
124 const std::vector<std::vector<int>>& it2iproj,
125 const std::vector<int>& iproj2l,
126 std::map<std::tuple<int, int, int, int>, int>& itiaiprojm2irow);
127 public:
136 RadialProjector(const bool realspace = false) {}
138
139 // it is more feasible to build interpolation table. this function will tabulate
140 // for those functions. This function will write the in-build tab_, to place the
141 // values of Jl[f](q) for each q and l.
142 // Here I provide two versions of tabulate, one for the case may be capable to
143 // avoid the memory copy operation, and the other one is for the case of
144 // std::vector<double> and std::vector<std::vector<double>>.
155 void _build_sbt_tab(const int nr,
156 const double* r,
157 const std::vector<double*>& radials,
158 const std::vector<int>& l,
159 const int nq, //< GlobalV::DQ
160 const double& dq); //< GlobalV::NQX
161 void _build_sbt_tab(const std::vector<double>& r,
162 const std::vector<std::vector<double>>& radials,
163 const std::vector<int>& l,
164 const int nq, //< GlobalV::DQ
165 const double& dq); //< GlobalV::NQX
166 // compatibility concern: for FS_Nonlocal_tools. Will not call sbtft so need omega
167 void _build_sbt_tab(const std::vector<int>& nproj,
168 const std::vector<double>& r,
169 const std::vector<std::vector<double>>& radials,
170 const std::vector<int>& l,
171 const int nq, //< GlobalV::DQ
172 const double& dq, //< GlobalV::NQX
173 const double& omega,
174 const int npol,
176 ModuleBase::matrix& nhtol);
189 void sbtft(const std::vector<ModuleBase::Vector3<double>>& qs,
190 std::vector<std::complex<double>>& out,
191 const char type = 'r', // 'r' for ket |>, 'l' for bra <|
192 const double& omega = 1.0,
193 const double& tpiba = 1.0); // 'n' for no gradient, 'x', 'y', 'z' for gradient in x, y, z direction
194
195 void sbfft(); // interface for SBFFT
196
197
198 private:
199 std::unique_ptr<ModuleBase::CubicSpline> cubspl_;
200 std::vector<int> l_;
201 };
202
233 void _mask_func(std::vector<double>& mask);
234
247 void _do_mask_on_radial(const int nr1,
248 const double* r,
249 const double* in,
250 const int nr2,
251 const double* mask,
252 double* out);
253}
254
255#endif // RADIAL_PROJECTION_H
3 elements vector
Definition vector3.h:22
Definition matrix.h:19
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:67
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:137
RadialProjector(const bool realspace=false)
Construct a new Radial Projector object.
Definition radial_proj.h:136
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:199
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:200
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:29
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