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