ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
cal_pLpR.h
Go to the documentation of this file.
1
86#include <cmath>
87#include <vector>
88#include <map>
89#include <tuple>
90#include <complex>
91#include <memory>
96
97namespace ModuleIO
98{
118 std::complex<double> cal_LzijR(
119 const std::unique_ptr<TwoCenterIntegrator>& calculator,
120 const int it, const int ia, const int il, const int iz, const int mi,
121 const int jt, const int ja, const int jl, const int jz, const int mj,
123
143 std::complex<double> cal_LyijR(
144 const std::unique_ptr<TwoCenterIntegrator>& calculator,
145 const int it, const int ia, const int il, const int iz, const int im,
146 const int jt, const int ja, const int jl, const int jz, const int jm,
148
168 std::complex<double> cal_LxijR(
169 const std::unique_ptr<TwoCenterIntegrator>& calculator,
170 const int it, const int ia, const int il, const int iz, const int im,
171 const int jt, const int ja, const int jl, const int jz, const int jm,
173
174 // the calculation of <phi_i|Lx/Ly/Lz|phi_j> matrix elements will be outputted
175 // in the way that indexed by:
176 // it, ia, il, iz, im, iRx, iRy, iRz, jt, ja, jl, jz, jm, in which the
177 // iRx, iRy, iRz are the indices of the supercell in which the two-center-integral
178 // it and jt are indexes of atomtypes,
179 // ia and ja are indexes of atoms within the atomtypes,
180 // il and jl are indexes of the angular momentum,
181 // iz and jz are indexes of the zeta functions
182 // im and jm are indexes of the magnetic quantum numbers.
183 // The output is a complex number, which is the value of the matrix element.
184 // Always the matrix is quite large, so direct print to file.
186 {
187 public:
188 // the default constructor is meaningless
203 const std::string& orbital_dir,
204 const UnitCell& ucell,
205 const double& search_radius,
206 const int tdestructor,
207 const int tgrid,
208 const int tatom,
209 const bool searchpbc,
210 std::ofstream* ptr_log = nullptr,
211 const int rank = 0);
213
214 void calculate(const std::string& prefix,
215 const std::string& outdir,
216 const UnitCell& ucell,
217 const int precision = 10,
218 const int rank = 0);
219
220 private:
221 // ofsrunning
222 std::ofstream* ofs_;
223 // the two-center-integrator
224 std::unique_ptr<TwoCenterIntegrator> calculator_;
225 // the spherical bessel transformer
227 // the radial collection
228 std::unique_ptr<RadialCollection> orb_;
229
230 // neighboring searcher
231 std::unique_ptr<Grid_Driver> neighbor_searcher_;
232
244 void kernel(std::ofstream* ofs,
245 const UnitCell& ucell,
246 const char dir = 'x',
247 const int precision = 10);
248 };
249} // namespace ModuleIO
A class that provides spherical Bessel transforms.
Definition spherical_bessel_transformer.h:69
3 elements vector
Definition vector3.h:24
Definition cal_pLpR.h:186
std::unique_ptr< TwoCenterIntegrator > calculator_
Definition cal_pLpR.h:224
ModuleBase::SphericalBesselTransformer sbt_
Definition cal_pLpR.h:226
std::unique_ptr< Grid_Driver > neighbor_searcher_
Definition cal_pLpR.h:231
void kernel(std::ofstream *ofs, const UnitCell &ucell, const char dir='x', const int precision=10)
calculate the <phi_i|Lx/Ly/Lz|phi_j> matrix elements. Due to the large size of the matrix,...
Definition cal_pLpR.cpp:256
std::ofstream * ofs_
Definition cal_pLpR.h:222
std::unique_ptr< RadialCollection > orb_
Definition cal_pLpR.h:228
Definition unitcell.h:15
Definition input_help.cpp:10
std::complex< double > cal_LyijR(const std::unique_ptr< TwoCenterIntegrator > &calculator, const int it, const int ia, const int il, const int iz, const int im, const int jt, const int ja, const int jl, const int jz, const int jm, const ModuleBase::Vector3< double > &vR)
calculate the <phi_i|Ly|phi_j> matrix elements, in which the Lz are the angular momentum operators,...
Definition cal_pLpR.cpp:117
std::complex< double > cal_LxijR(const std::unique_ptr< TwoCenterIntegrator > &calculator, const int it, const int ia, const int il, const int iz, const int im, const int jt, const int ja, const int jl, const int jz, const int jm, const ModuleBase::Vector3< double > &vR)
calculate the <phi_i|Lx|phi_j> matrix elements, in which the Lz are the angular momentum operators,...
Definition cal_pLpR.cpp:63
std::complex< double > cal_LzijR(const std::unique_ptr< TwoCenterIntegrator > &calculator, const int it, const int ia, const int il, const int iz, const int mi, const int jt, const int ja, const int jl, const int jz, const int mj, const ModuleBase::Vector3< double > &vR)
calculate the <phi_i|Lz|phi_j> matrix elements, in which the Lz are the angular momentum operators,...
Definition cal_pLpR.cpp:49
void calculate()
Definition main.cpp:24