ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
dftu.h
Go to the documentation of this file.
1#ifndef DFTU_H
2#define DFTU_H
3
4#include "source_cell/klist.h"
7#ifdef __LCAO
13#include "source_lcao/force_stress_arrays.h" // mohan add 2024-06-15
14#endif
15
16#include <string>
17#include <vector>
18
19//==========================================================
20// CLASS :
21// NAME : DTFU (DFT+U)
22//==========================================================
23namespace ModuleDFTU
24{
25
26class DFTU
27{
28
29 public:
30 DFTU(); // constructor
31 ~DFTU(); // deconstructor
32
33 //=============================================================
34 // In dftu.cpp
35 // Initialization & Calculating energy
36 //=============================================================
37 public:
38 // allocate relevant data strcutures
39 void init(UnitCell& cell, // unitcell class
40 const Parallel_Orbitals* pv,
41 const int nks
42#ifdef __LCAO
43 , const LCAO_Orbitals* orb = nullptr
44#endif
45 );
46
47 static DFTU* get_instance();
48
49 // calculate the energy correction
50 void cal_energy_correction(const UnitCell& ucell, const int istep);
51 double get_energy(){return EU;}
52 void uramping_update(); // update U by uramping
53 bool u_converged(); // check if U is converged
54
55 std::vector<double> U = {}; // U (Hubbard parameter U)
56 std::vector<double> U0; // U0 (target Hubbard parameter U0)
57 std::vector<int> orbital_corr = {}; //
58 double uramping; // increase U by uramping, default is -1.0
59 int omc; // occupation matrix control
60 int mixing_dftu; //whether to mix locale
61
62 double EU; //+U energy
63 private:
64 const Parallel_Orbitals* paraV = nullptr;
65 int cal_type = 3; // 1:dftu_tpye=1, dc=1; 2:dftu_type=1, dc=2; 3:dftu_tpye=2, dc=1; 4:dftu_tpye=2, dc=2;
66
67 // FIXME: the following variable does not have static lifetime;
68 // while the present class is used via a global variable. This has
69 // potential to cause dangling pointer issues.
70#ifdef __LCAO
71 const LCAO_Orbitals* ptr_orb_ = nullptr;
72 std::vector<double> orb_cutoff_;
73#endif
74
75 // transform between iwt index and it, ia, L, N and m index
76 std::vector<std::vector<std::vector<std::vector<std::vector<int>>>>>
77 iatlnmipol2iwt; // iatlnm2iwt[iat][l][n][m][ipol]
78
79#ifdef __LCAO
80 //=============================================================
81 // In dftu_hamilt.cpp
82 // For calculating contribution to Hamiltonian matrices
83 //=============================================================
84 public:
85 void cal_eff_pot_mat_complex(const int ik, std::complex<double>* eff_pot, const std::vector<int>& isk, const std::complex<double>* sk);
86 void cal_eff_pot_mat_real(const int ik, double* eff_pot, const std::vector<int>& isk, const double* sk);
87 void cal_eff_pot_mat_R_double(const int ispin, double* SR, double* HR);
88 void cal_eff_pot_mat_R_complex_double(const int ispin, std::complex<double>* SR, std::complex<double>* HR);
89#endif
90
91 //=============================================================
92 // In dftu_occup.cpp
93 // For calculating occupation matrix and saving to locale
94 // and other operations of locale: copy,zero out,mix
95 //=============================================================
96 public:
99 void cal_occ_pw(const int iter, const void* psi_in, const ModuleBase::matrix& wg_in, const UnitCell& cell, const double& mixing_beta);
101 void cal_VU_pot_pw(const int spin);
103 const std::complex<double>* get_eff_pot_pw(const int iat) const { return &(eff_pot_pw[this->eff_pot_pw_index[iat]]); }
104 int get_size_eff_pot_pw() const { return eff_pot_pw.size(); }
105
106#ifdef __LCAO
107 // calculate the local occupation number matrix
108 void cal_occup_m_k(const int iter,
109 const UnitCell& ucell,
110 const std::vector<std::vector<std::complex<double>>>& dm_k,
111 const K_Vectors& kv,
112 const double& mixing_beta,
113 hamilt::Hamilt<std::complex<double>>* p_ham);
114 void cal_occup_m_gamma(const int iter,
115 const UnitCell& ucell,
116 const std::vector<std::vector<double>>& dm_gamma,
117 const double& mixing_beta,
119#endif
120
121 // dftu can be calculated only after locale has been initialed
122 bool initialed_locale = false;
123
124 private:
125 void copy_locale(const UnitCell& ucell);
126 void zero_locale(const UnitCell& ucell);
127 void mix_locale(const UnitCell& ucell,const double& mixing_beta);
128
129 std::vector<std::complex<double>> eff_pot_pw;
130 std::vector<int> eff_pot_pw_index;
131
132public:
133 // local occupancy matrix of the correlated subspace
134 // locale: the out put local occupation number matrix of correlated electrons in the current electronic step
135 // locale_save: the input local occupation number matrix of correlated electrons in the current electronic step
136 std::vector<std::vector<std::vector<std::vector<ModuleBase::matrix>>>> locale; // locale[iat][l][n][spin](m1,m2)
137 std::vector<std::vector<std::vector<std::vector<ModuleBase::matrix>>>> locale_save; // locale_save[iat][l][n][spin](m1,m2)
138#ifdef __LCAO
139private:
140 //=============================================================
141 // In dftu_tools.cpp
142 // For calculating onsite potential, which is used
143 // for both Hamiltonian and force/stress
144 //=============================================================
145
146 void cal_VU_pot_mat_complex(const int spin, const bool newlocale, std::complex<double>* VU);
147 void cal_VU_pot_mat_real(const int spin, const bool newlocale, double* VU);
148
149 double get_onebody_eff_pot(const int T,
150 const int iat,
151 const int L,
152 const int N,
153 const int spin,
154 const int m0,
155 const int m1,
156 const bool newlocale);
157
158 //=============================================================
159 // In dftu_folding.cpp
160 // Subroutines for folding S and dS matrix
161 //=============================================================
162
163 void fold_dSR_gamma(const UnitCell& ucell,
164 const Parallel_Orbitals& pv,
165 const Grid_Driver* gd,
166 double* dsloc_x,
167 double* dsloc_y,
168 double* dsloc_z,
169 double* dh_r,
170 const int dim1,
171 const int dim2,
172 double* dSR_gamma);
173
174 // dim = 0 : S, for Hamiltonian
175 // dim = 1-3 : dS, for force
176 // dim = 4-6 : dS * dR, for stress
177
178 void folding_matrix_k(const UnitCell& ucell,
179 const Grid_Driver& gd,
181 const Parallel_Orbitals& pv,
182 const int ik,
183 const int dim1,
184 const int dim2,
185 std::complex<double>* mat_k,
186 const ModuleBase::Vector3<double>& kvec_d);
187
193 void folding_matrix_k_new(const int ik,
194 hamilt::Hamilt<std::complex<double>>* p_ham);
195
196 //=============================================================
197 // In dftu_force.cpp
198 // For calculating force and stress fomr DFT+U
199 //=============================================================
200 public:
201 void force_stress(const UnitCell& ucell,
202 const Grid_Driver& gd,
203 const elecstate::ElecState* pelec,
204 const Parallel_Orbitals& pv,
206 ModuleBase::matrix& force_dftu,
207 ModuleBase::matrix& stress_dftu,
208 const K_Vectors& kv);
209
210 private:
211 void cal_force_k(const UnitCell& ucell,
212 const Grid_Driver& gd,
214 const Parallel_Orbitals& pv,
215 const int ik,
216 const std::complex<double>* rho_VU,
217 ModuleBase::matrix& force_dftu,
218 const ModuleBase::Vector3<double>& kvec_d);
219
220 void cal_stress_k(const UnitCell& ucell,
221 const Grid_Driver& gd,
223 const Parallel_Orbitals& pv,
224 const int ik,
225 const std::complex<double>* rho_VU,
226 ModuleBase::matrix& stress_dftu,
227 const ModuleBase::Vector3<double>& kvec_d);
228
229 void cal_force_gamma(const UnitCell& ucell,
230 const double* rho_VU,
231 const Parallel_Orbitals& pv,
232 double* dsloc_x,
233 double* dsloc_y,
234 double* dsloc_z,
235 ModuleBase::matrix& force_dftu);
236
237 void cal_stress_gamma(const UnitCell& ucell,
238 const Parallel_Orbitals& pv,
239 const Grid_Driver* gd,
240 double* dsloc_x,
241 double* dsloc_y,
242 double* dsloc_z,
243 double* dh_r,
244 const double* rho_VU,
245 ModuleBase::matrix& stress_dftu);
246#endif
247
248 //=============================================================
249 // In dftu_io.cpp
250 // For reading/writing/broadcasting/copying relevant data structures
251 //=============================================================
252 public:
253 void output(const UnitCell& ucell);
254
255 private:
256 void write_occup_m(const UnitCell& ucell,
257 std::ofstream& ofs,
258 bool diag=false);
259 void read_occup_m(const UnitCell& ucell,
260 const std::string& fn);
261 void local_occup_bcast(const UnitCell& ucell);
262
263 //=============================================================
264 // In dftu_yukawa.cpp
265 // Relevant for calculating U using Yukawa potential
266 //=============================================================
267
268 public:
269 bool Yukawa; // 1:use Yukawa potential; 0: do not use Yukawa potential
270 void cal_slater_UJ(const UnitCell& ucell, double** rho, const int& nrxx);
271
272 private:
273 double lambda; // the parameter in Yukawa potential
274 std::vector<std::vector<std::vector<std::vector<double>>>> Fk; // slater integral:Fk[T][L][N][k]
275 std::vector<std::vector<std::vector<double>>> U_Yukawa; // U_Yukawa[T][L][N]
276 std::vector<std::vector<std::vector<double>>> J_Yukawa; // J_Yukawa[T][L][N]
277
278 void cal_slater_Fk(const UnitCell& ucell,const int L, const int T); // L:angular momnet, T:atom type
279 void cal_yukawa_lambda(double** rho, const int& nrxx);
280
281 double spherical_Bessel(const int k, const double r, const double lambda);
282 double spherical_Hankel(const int k, const double r, const double lambda);
283
284#ifdef __LCAO
285 public:
291 const hamilt::HContainer<double>* get_dmr(int ispin) const;
296 void set_dmr(const elecstate::DensityMatrix<double, double>* dm_in_dftu_d);
297 void set_dmr(const elecstate::DensityMatrix<std::complex<double>, double>* dm_in_dftu_cd);
298
299 private:
300 const UnitCell* ucell = nullptr;
301 const elecstate::DensityMatrix<double, double>* dm_in_dftu_d = nullptr;
302 const elecstate::DensityMatrix<std::complex<double>, double>* dm_in_dftu_cd = nullptr;
303#endif
304};
305
306#ifdef __LCAO
307template <typename T>
308void dftu_cal_occup_m(const int iter,
309 const UnitCell& ucell,
310 const std::vector<std::vector<T>>& dm,
311 const K_Vectors& kv,
312 const double& mixing_beta,
313 hamilt::Hamilt<T>* p_ham);
314#endif
315
316} // namespace ModuleDFTU
317
318namespace GlobalC
319{
320 extern ModuleDFTU::DFTU dftu;
321}
322#endif
Definition force_stress_arrays.h:5
Definition sltk_grid_driver.h:43
Definition klist.h:13
Definition ORB_read.h:19
3 elements vector
Definition vector3.h:22
Definition matrix.h:19
Definition dftu.h:27
void uramping_update()
Definition dftu.cpp:383
void local_occup_bcast(const UnitCell &ucell)
Definition dftu_io.cpp:390
int omc
Definition dftu.h:59
const Parallel_Orbitals * paraV
Definition dftu.h:64
bool Yukawa
Definition dftu.h:269
int get_size_eff_pot_pw() const
Definition dftu.h:104
void zero_locale(const UnitCell &ucell)
Definition dftu_occup.cpp:50
std::vector< std::vector< std::vector< std::vector< ModuleBase::matrix > > > > locale_save
Definition dftu.h:137
void cal_slater_UJ(const UnitCell &ucell, double **rho, const int &nrxx)
bool initialed_locale
Definition dftu.h:122
void write_occup_m(const UnitCell &ucell, std::ofstream &ofs, bool diag=false)
Definition dftu_io.cpp:76
DFTU()
Definition for_testing_input_conv.h:139
std::vector< std::vector< std::vector< std::vector< ModuleBase::matrix > > > > locale
Definition dftu.h:136
~DFTU()
Definition for_testing_input_conv.h:142
std::vector< std::vector< std::vector< std::vector< std::vector< int > > > > > iatlnmipol2iwt
Definition dftu.h:77
std::vector< std::vector< std::vector< std::vector< double > > > > Fk
Definition dftu.h:274
static DFTU * get_instance()
Definition dftu_pw.cpp:10
void read_occup_m(const UnitCell &ucell, const std::string &fn)
Definition dftu_io.cpp:225
void init(UnitCell &cell, const Parallel_Orbitals *pv, const int nks)
Definition dftu.cpp:38
double uramping
Definition dftu.h:58
int mixing_dftu
Definition dftu.h:60
double spherical_Bessel(const int k, const double r, const double lambda)
std::vector< double > U0
Definition dftu.h:56
void cal_slater_Fk(const UnitCell &ucell, const int L, const int T)
int cal_type
Definition dftu.h:65
double get_energy()
Definition dftu.h:51
std::vector< std::complex< double > > eff_pot_pw
Definition dftu.h:129
void copy_locale(const UnitCell &ucell)
Definition dftu_occup.cpp:12
std::vector< double > U
Definition dftu.h:55
double lambda
Definition dftu.h:273
std::vector< int > eff_pot_pw_index
Definition dftu.h:130
void cal_VU_pot_pw(const int spin)
calculate the local DFT+U effective potential matrix for PW base.
Definition dftu_pw.cpp:207
std::vector< std::vector< std::vector< double > > > J_Yukawa
Definition dftu.h:276
void cal_occ_pw(const int iter, const void *psi_in, const ModuleBase::matrix &wg_in, const UnitCell &cell, const double &mixing_beta)
calculate occupation matrix for DFT+U
Definition dftu_pw.cpp:15
void mix_locale(const UnitCell &ucell, const double &mixing_beta)
Definition dftu_occup.cpp:88
void cal_yukawa_lambda(double **rho, const int &nrxx)
std::vector< int > orbital_corr
Definition dftu.h:57
double spherical_Hankel(const int k, const double r, const double lambda)
const std::complex< double > * get_eff_pot_pw(const int iat) const
get effective potential matrix for PW base
Definition dftu.h:103
void cal_energy_correction(const UnitCell &ucell, const int istep)
bool u_converged()
Definition dftu.cpp:403
std::vector< std::vector< std::vector< double > > > U_Yukawa
Definition dftu.h:275
double EU
Definition dftu.h:62
Definition parallel_orbitals.h:9
Definition unitcell.h:16
Definition density_matrix.h:36
Definition elecstate.h:15
Definition hcontainer.h:144
Definition hamilt.h:16
Definition output.h:13
#define N
Definition exp.cpp:24
#define T
Definition exp.cpp:237
Definition cal_epsilon_test.cpp:31
ModuleDFTU::DFTU dftu
Definition for_testing_input_conv.h:254
Definition dftu.cpp:29
base device SOURCES math_dngvd_test cpp endif() if(ENABLE_GOOGLEBENCH) AddTest(TARGET PERF_MODULE_HSOLVER_KERNELS LIBS parameter $
Definition CMakeLists.txt:10