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
61#include <cmath>
62#include <vector>
63#include <map>
64#include <tuple>
65#include <complex>
66#include <memory>
71
72namespace ModuleIO
73{
93 std::complex<double> cal_LzijR(
94 const std::unique_ptr<TwoCenterIntegrator>& calculator,
95 const int it, const int ia, const int il, const int iz, const int mi,
96 const int jt, const int ja, const int jl, const int jz, const int mj,
98
118 std::complex<double> cal_LyijR(
119 const std::unique_ptr<TwoCenterIntegrator>& calculator,
120 const int it, const int ia, const int il, const int iz, const int im,
121 const int jt, const int ja, const int jl, const int jz, const int jm,
123
143 std::complex<double> cal_LxijR(
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
149 // the calculation of <phi_i|Lx/Ly/Lz|phi_j> matrix elements will be outputted
150 // in the way that indexed by:
151 // it, ia, il, iz, im, iRx, iRy, iRz, jt, ja, jl, jz, jm, in which the
152 // iRx, iRy, iRz are the indices of the supercell in which the two-center-integral
153 // it and jt are indexes of atomtypes,
154 // ia and ja are indexes of atoms within the atomtypes,
155 // il and jl are indexes of the angular momentum,
156 // iz and jz are indexes of the zeta functions
157 // im and jm are indexes of the magnetic quantum numbers.
158 // The output is a complex number, which is the value of the matrix element.
159 // Always the matrix is quite large, so direct print to file.
161 {
162 public:
163 // the default constructor is meaningless
178 const std::string& orbital_dir,
179 const UnitCell& ucell,
180 const double& search_radius,
181 const int tdestructor,
182 const int tgrid,
183 const int tatom,
184 const bool searchpbc,
185 std::ofstream* ptr_log = nullptr,
186 const int rank = 0);
188
189 void calculate(const std::string& prefix,
190 const std::string& outdir,
191 const UnitCell& ucell,
192 const int precision = 10,
193 const int rank = 0);
194
195 private:
196 // ofsrunning
197 std::ofstream* ofs_;
198 // the two-center-integrator
199 std::unique_ptr<TwoCenterIntegrator> calculator_;
200 // the spherical bessel transformer
202 // the radial collection
203 std::unique_ptr<RadialCollection> orb_;
204
205 // neighboring searcher
206 std::unique_ptr<Grid_Driver> neighbor_searcher_;
207
219 void kernel(std::ofstream* ofs,
220 const UnitCell& ucell,
221 const char dir = 'x',
222 const int precision = 10);
223 };
224} // namespace ModuleIO
A class that provides spherical Bessel transforms.
Definition spherical_bessel_transformer.h:69
3 elements vector
Definition vector3.h:22
Definition cal_pLpR.h:161
std::unique_ptr< TwoCenterIntegrator > calculator_
Definition cal_pLpR.h:199
ModuleBase::SphericalBesselTransformer sbt_
Definition cal_pLpR.h:201
std::unique_ptr< Grid_Driver > neighbor_searcher_
Definition cal_pLpR.h:206
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:168
std::ofstream * ofs_
Definition cal_pLpR.h:197
std::unique_ptr< RadialCollection > orb_
Definition cal_pLpR.h:203
Definition unitcell.h:16
This class has two functions: restart psi from the previous calculation, and write psi to the disk.
Definition cal_dos.h:9
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:46
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:69
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:35
void calculate()
Definition main.cpp:24