ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
Exx_LRI.hpp
Go to the documentation of this file.
1//=======================
2// AUTHOR : Peize Lin
4// DATE : 2022-08-17
5//=======================
6
7#ifndef EXX_LRI_HPP
8#define EXX_LRI_HPP
9
10#include "Exx_LRI.h"
11#include "RI_2D_Comm.h"
12#include "RI_Util.h"
17#include "source_base/timer.h"
21
22#include <RI/distribute/Distribute_Equally.h>
23#include <RI/global/Map_Operator-3.h>
24
25#include <fstream>
26#include <string>
27
28template<typename Tdata>
29void Exx_LRI<Tdata>::init(const MPI_Comm &mpi_comm_in,
30 const UnitCell &ucell,
31 const K_Vectors &kv_in,
32 const LCAO_Orbitals& orb)
33{
34 ModuleBase::TITLE("Exx_LRI","init");
35 ModuleBase::timer::tick("Exx_LRI", "init");
36
37 this->mpi_comm = mpi_comm_in;
38 this->p_kv = &kv_in;
39 this->orb_cutoff_ = orb.cutoffs();
40
41 this->lcaos = Exx_Abfs::Construct_Orbs::change_orbs( orb, this->info.kmesh_times );
42
43 const std::vector<std::vector<std::vector<Numerical_Orbital_Lm>>>
44 abfs_same_atom = Exx_Abfs::Construct_Orbs::abfs_same_atom(ucell, orb, this->lcaos, this->info.kmesh_times, this->info.pca_threshold );
45 if(this->info.files_abfs.empty())
46 { this->abfs = abfs_same_atom;}
47 else
48 { this->abfs = Exx_Abfs::IO::construct_abfs( abfs_same_atom, orb, this->info.files_abfs, this->info.kmesh_times ); }
50
51 for( size_t T=0; T!=this->abfs.size(); ++T )
52 { GlobalC::exx_info.info_ri.abfs_Lmax = std::max( GlobalC::exx_info.info_ri.abfs_Lmax, static_cast<int>(this->abfs[T].size())-1 ); }
53
54 this->coulomb_settings = RI_Util::update_coulomb_settings(this->info.coulomb_param, ucell, this->p_kv);
55
56 bool init_MGT = true;
57 for(const auto &settings_list : this->coulomb_settings)
58 {
59 this->exx_objs[settings_list.first].abfs_ccp = Conv_Coulomb_Pot_K::cal_orbs_ccp(this->abfs, settings_list.second.second, this->info.ccp_rmesh_times);
60 this->exx_objs[settings_list.first].cv.set_orbitals(ucell, orb,
61 this->lcaos, this->abfs, this->exx_objs[settings_list.first].abfs_ccp,
62 this->info.kmesh_times, this->MGT, init_MGT, settings_list.second.first );
63 init_MGT = false; // only init once
64 }
65
66 ModuleBase::timer::tick("Exx_LRI", "init");
67}
68
69template<typename Tdata>
71 const bool write_cv)
72{
73 ModuleBase::TITLE("Exx_LRI","cal_exx_ions");
74 ModuleBase::timer::tick("Exx_LRI", "cal_exx_ions");
75
76 // init_radial_table_ions( cal_atom_centres_core(atom_pairs_core_origin), atom_pairs_core_origin );
77
78 // this->m_abfsabfs.init_radial_table(Rradial);
79 // this->m_abfslcaos_lcaos.init_radial_table(Rradial);
80
81 std::vector<TA> atoms(ucell.nat);
82 for(int iat=0; iat<ucell.nat; ++iat)
83 { atoms[iat] = iat; }
84 std::map<TA,TatomR> atoms_pos;
85 for(int iat=0; iat<ucell.nat; ++iat)
86 { atoms_pos[iat] = RI_Util::Vector3_to_array3( ucell.atoms[ ucell.iat2it[iat] ].tau[ ucell.iat2ia[iat] ] ); }
87 const std::array<TatomR,Ndim> latvec
91 const std::array<Tcell,Ndim> period = {this->p_kv->nmp[0], this->p_kv->nmp[1], this->p_kv->nmp[2]};
92
93 this->exx_lri.set_parallel(this->mpi_comm, atoms_pos, latvec, period);
94
95 // std::max(3) for gamma_only, list_A2 should contain cell {-1,0,1}. In the future distribute will be neighbour.
96 const std::array<Tcell,Ndim> period_Vs = LRI_CV_Tools::cal_latvec_range<Tcell>(1+this->info.ccp_rmesh_times, ucell, orb_cutoff_);
97 const std::pair<std::vector<TA>, std::vector<std::vector<std::pair<TA,std::array<Tcell,Ndim>>>>>
98 list_As_Vs = RI::Distribute_Equally::distribute_atoms_periods(this->mpi_comm, atoms, period_Vs, 2, false);
99
100 std::map<TA,std::map<TAC,RI::Tensor<Tdata>>> Vs;
101 std::map<TA, std::map<TAC, std::array<RI::Tensor<Tdata>, Ndim>>> dVs;
102 for(const auto &settings_list : this->coulomb_settings)
103 {
104 std::map<TA,std::map<TAC,RI::Tensor<Tdata>>>
105 Vs_temp = this->exx_objs[settings_list.first].cv.cal_Vs(ucell,
106 list_As_Vs.first, list_As_Vs.second[0],
107 {{"writable_Vws",true}});
108 this->exx_objs[settings_list.first].cv.Vws = LRI_CV_Tools::get_CVws(ucell,Vs_temp);
109 Vs = Vs.empty() ? Vs_temp : LRI_CV_Tools::add(Vs, Vs_temp);
110
112 {
113 std::map<TA, std::map<TAC, std::array<RI::Tensor<Tdata>, Ndim>>>
114 dVs_temp = this->exx_objs[settings_list.first].cv.cal_dVs(ucell,
115 list_As_Vs.first, list_As_Vs.second[0],
116 {{"writable_dVws",true}});
117 this->exx_objs[settings_list.first].cv.dVws = LRI_CV_Tools::get_dCVws(ucell,dVs_temp);
118 dVs = dVs.empty() ? dVs_temp : LRI_CV_Tools::add(dVs, dVs_temp);
119 }
120 }
121 if (write_cv && GlobalV::MY_RANK == 0)
123 this->exx_lri.set_Vs(std::move(Vs), this->info.V_threshold);
124
126 {
127 std::array<std::map<TA, std::map<TAC, RI::Tensor<Tdata>>>, Ndim>
128 dVs_order = LRI_CV_Tools::change_order(std::move(dVs));
129 this->exx_lri.set_dVs(std::move(dVs_order), this->info.V_grad_threshold);
131 {
132 std::array<std::array<std::map<TA,std::map<TAC,RI::Tensor<Tdata>>>,3>,3> dVRs = LRI_CV_Tools::cal_dMRs(ucell,dVs_order);
133 this->exx_lri.set_dVRs(std::move(dVRs), this->info.V_grad_R_threshold);
134 }
135 }
136
137 const std::array<Tcell,Ndim> period_Cs = LRI_CV_Tools::cal_latvec_range<Tcell>(2, ucell,orb_cutoff_);
138 const std::pair<std::vector<TA>, std::vector<std::vector<std::pair<TA,std::array<Tcell,Ndim>>>>>
139 list_As_Cs = RI::Distribute_Equally::distribute_atoms_periods(this->mpi_comm, atoms, period_Cs, 2, false);
140
141 std::map<TA,std::map<TAC,RI::Tensor<Tdata>>> Cs;
142 std::map<TA, std::map<TAC, std::array<RI::Tensor<Tdata>, 3>>> dCs;
143 for(const auto &settings_list : this->coulomb_settings)
144 {
145 if(settings_list.second.first)
146 {
147 std::pair<std::map<TA, std::map<TAC, RI::Tensor<Tdata>>>,
148 std::map<TA, std::map<TAC, std::array<RI::Tensor<Tdata>, 3>>>>
149 Cs_dCs = this->exx_objs[settings_list.first].cv.cal_Cs_dCs(
150 ucell,
151 list_As_Cs.first, list_As_Cs.second[0],
152 {{"cal_dC",PARAM.inp.cal_force||PARAM.inp.cal_stress},
153 {"writable_Cws",true}, {"writable_dCws",true}, {"writable_Vws",false}, {"writable_dVws",false}});
154 std::map<TA,std::map<TAC,RI::Tensor<Tdata>>> &Cs_temp = std::get<0>(Cs_dCs);
155 this->exx_objs[settings_list.first].cv.Cws = LRI_CV_Tools::get_CVws(ucell,Cs_temp);
156 Cs = Cs.empty() ? Cs_temp : LRI_CV_Tools::add(Cs, Cs_temp);
157
159 {
160 std::map<TA, std::map<TAC, std::array<RI::Tensor<Tdata>, 3>>> &dCs_temp = std::get<1>(Cs_dCs);
161 this->exx_objs[settings_list.first].cv.dCws = LRI_CV_Tools::get_dCVws(ucell,dCs_temp);
162 dCs = dCs.empty() ? dCs_temp : LRI_CV_Tools::add(dCs, dCs_temp);
163 }
164 }
165 }
166 if (write_cv && GlobalV::MY_RANK == 0)
168 this->exx_lri.set_Cs(std::move(Cs), this->info.C_threshold);
169
171 {
172 std::array<std::map<TA, std::map<TAC, RI::Tensor<Tdata>>>, Ndim>
173 dCs_order = LRI_CV_Tools::change_order(std::move(dCs));
174 this->exx_lri.set_dCs(std::move(dCs_order), this->info.C_grad_threshold);
176 {
177 std::array<std::array<std::map<TA,std::map<TAC,RI::Tensor<Tdata>>>,3>,3> dCRs = LRI_CV_Tools::cal_dMRs(ucell,dCs_order);
178 this->exx_lri.set_dCRs(std::move(dCRs), this->info.C_grad_R_threshold);
179 }
180 }
181 ModuleBase::timer::tick("Exx_LRI", "cal_exx_ions");
182}
183
184template<typename Tdata>
185void Exx_LRI<Tdata>::cal_exx_elec(const std::vector<std::map<TA, std::map<TAC, RI::Tensor<Tdata>>>>& Ds,
186 const UnitCell& ucell,
187 const Parallel_Orbitals& pv,
188 const ModuleSymmetry::Symmetry_rotation* p_symrot)
189{
190 ModuleBase::TITLE("Exx_LRI","cal_exx_elec");
191 ModuleBase::timer::tick("Exx_LRI", "cal_exx_elec");
192
193 const std::vector<std::tuple<std::set<TA>, std::set<TA>>> judge = RI_2D_Comm::get_2D_judge(ucell,pv);
194
195 if(p_symrot)
196 { this->exx_lri.set_symmetry(true, p_symrot->get_irreducible_sector()); }
197 else
198 { this->exx_lri.set_symmetry(false, {}); }
199
200 this->Hexxs.resize(PARAM.inp.nspin);
201 this->Eexx = 0;
202 for(int is=0; is<PARAM.inp.nspin; ++is)
203 {
204 const std::string suffix = ((PARAM.inp.cal_force || PARAM.inp.cal_stress) ? std::to_string(is) : "");
205
206 this->exx_lri.set_Ds(Ds[is], this->info.dm_threshold, suffix);
207 this->exx_lri.cal_Hs({ "","",suffix });
208
209 if (!p_symrot)
210 {
211 this->Hexxs[is] = RI::Communicate_Tensors_Map_Judge::comm_map2_first(
212 this->mpi_comm, std::move(this->exx_lri.Hs), std::get<0>(judge[is]), std::get<1>(judge[is]));
213 }
214 else
215 {
216 // reduce but not repeat
217 auto Hs_a2D = this->exx_lri.post_2D.set_tensors_map2(this->exx_lri.Hs);
218 // rotate locally without repeat
219 Hs_a2D = p_symrot->restore_HR(ucell.symm, ucell.atoms, ucell.st, 'H', Hs_a2D);
220 // cal energy using full Hs without repeat
221 this->exx_lri.energy = this->exx_lri.post_2D.cal_energy(
222 this->exx_lri.post_2D.saves["Ds_" + suffix],
223 this->exx_lri.post_2D.set_tensors_map2(Hs_a2D));
224 // get repeated full Hs for abacus
225 this->Hexxs[is] = RI::Communicate_Tensors_Map_Judge::comm_map2_first(
226 this->mpi_comm, std::move(Hs_a2D), std::get<0>(judge[is]), std::get<1>(judge[is]));
227 }
228 this->Eexx += std::real(this->exx_lri.energy);
229 post_process_Hexx(this->Hexxs[is]);
230 }
231 this->Eexx = post_process_Eexx(this->Eexx);
232 this->exx_lri.set_symmetry(false, {});
233 ModuleBase::timer::tick("Exx_LRI", "cal_exx_elec");
234}
235
236template<typename Tdata>
237void Exx_LRI<Tdata>::post_process_Hexx( std::map<TA, std::map<TAC, RI::Tensor<Tdata>>> &Hexxs_io ) const
238{
239 ModuleBase::TITLE("Exx_LRI","post_process_Hexx");
240 constexpr Tdata frac = -1 * 2; // why? Hartree to Ry?
241 const std::function<void(RI::Tensor<Tdata>&)>
242 multiply_frac = [&frac](RI::Tensor<Tdata> &t)
243 { t = t*frac; };
244 RI::Map_Operator::for_each( Hexxs_io, multiply_frac );
245}
246
247template<typename Tdata>
248double Exx_LRI<Tdata>::post_process_Eexx(const double& Eexx_in) const
249{
250 ModuleBase::TITLE("Exx_LRI","post_process_Eexx");
251 const double SPIN_multiple = std::map<int, double>{ {1,2}, {2,1}, {4,1} }.at(PARAM.inp.nspin); // why?
252 const double frac = -SPIN_multiple;
253 return frac * Eexx_in;
254}
255
256/*
257post_process_old
258{
259 // D
260 const std::map<int,double> SPIN_multiple = {{1,0.5}, {2,1}, {4,1}}; // ???
261 DR *= SPIN_multiple.at(NSPIN);
262
263 // H
264 HR *= -2;
265
266 // E
267 const std::map<int,double> SPIN_multiple = {{1,2}, {2,1}, {4,1}}; // ???
268 energy *= SPIN_multiple.at(PARAM.inp.nspin); // ?
269 energy /= 2; // /2 for Ry
270}
271*/
272
273template<typename Tdata>
275{
276 ModuleBase::TITLE("Exx_LRI","cal_exx_force");
277 ModuleBase::timer::tick("Exx_LRI", "cal_exx_force");
278
279 this->force_exx.create(nat, Ndim);
280 for(int is=0; is<PARAM.inp.nspin; ++is)
281 {
282 this->exx_lri.cal_force({"","",std::to_string(is),"",""});
283 for(std::size_t idim=0; idim<Ndim; ++idim) {
284 for(const auto &force_item : this->exx_lri.force[idim]) {
285 this->force_exx(force_item.first, idim) += std::real(force_item.second);
286 } }
287 }
288
289 const double SPIN_multiple = std::map<int,double>{{1,2}, {2,1}, {4,1}}.at(PARAM.inp.nspin); // why?
290 const double frac = -2 * SPIN_multiple; // why?
291 this->force_exx *= frac;
292 ModuleBase::timer::tick("Exx_LRI", "cal_exx_force");
293}
294
295
296template<typename Tdata>
297void Exx_LRI<Tdata>::cal_exx_stress(const double& omega, const double& lat0)
298{
299 ModuleBase::TITLE("Exx_LRI","cal_exx_stress");
300 ModuleBase::timer::tick("Exx_LRI", "cal_exx_stress");
301
302 this->stress_exx.create(Ndim, Ndim);
303 for(int is=0; is<PARAM.inp.nspin; ++is)
304 {
305 this->exx_lri.cal_stress({"","",std::to_string(is),"",""});
306 for(std::size_t idim0=0; idim0<Ndim; ++idim0) {
307 for(std::size_t idim1=0; idim1<Ndim; ++idim1) {
308 this->stress_exx(idim0,idim1) += std::real(this->exx_lri.stress(idim0,idim1));
309 } }
310 }
311
312 const double SPIN_multiple = std::map<int,double>{{1,2}, {2,1}, {4,1}}.at(PARAM.inp.nspin); // why?
313 const double frac = 2 * SPIN_multiple / omega * lat0; // why?
314 this->stress_exx *= frac;
315
316 ModuleBase::timer::tick("Exx_LRI", "cal_exx_stress");
317}
318
319/*
320template<typename Tdata>
321std::vector<std::vector<int>> Exx_LRI<Tdata>::get_abfs_nchis() const
322{
323 std::vector<std::vector<int>> abfs_nchis;
324 for (const auto& abfs_T : this->abfs)
325 {
326 std::vector<int> abfs_nchi_T;
327 for (const auto& abfs_L : abfs_T)
328 { abfs_nchi_T.push_back(abfs_L.size()); }
329 abfs_nchis.push_back(abfs_nchi_T);
330 }
331 return abfs_nchis;
332}
333*/
334
335#endif
std::vector< ModuleBase::Vector3< double > > tau
Definition atom_spec.h:36
static void print_orbs_size(const UnitCell &ucell, const std::vector< std::vector< std::vector< Numerical_Orbital_Lm > > > &orbs, std::ostream &os)
Definition exx_abfs-construct_orbs.cpp:474
static std::vector< std::vector< std::vector< Numerical_Orbital_Lm > > > abfs_same_atom(const UnitCell &ucell, const LCAO_Orbitals &orb, const std::vector< std::vector< std::vector< Numerical_Orbital_Lm > > > &lcaos, const double kmesh_times_mot, const double times_threshold=0)
Definition exx_abfs-construct_orbs.cpp:81
static std::vector< std::vector< std::vector< Numerical_Orbital_Lm > > > change_orbs(const LCAO_Orbitals &orb_in, const double kmesh_times)
Definition exx_abfs-construct_orbs.cpp:11
static std::vector< std::vector< std::vector< Numerical_Orbital_Lm > > > construct_abfs(const LCAO_Orbitals &orbs, const std::vector< std::string > &files_abfs, const double kmesh_times=1)
Definition exx_abfs-io.cpp:14
void post_process_Hexx(std::map< TA, std::map< TAC, RI::Tensor< Tdata > > > &Hexxs_io) const
Definition Exx_LRI.hpp:237
void init(const MPI_Comm &mpi_comm_in, const UnitCell &ucell, const K_Vectors &kv_in, const LCAO_Orbitals &orb)
Definition Exx_LRI.hpp:29
void cal_exx_elec(const std::vector< std::map< TA, std::map< TAC, RI::Tensor< Tdata > > > > &Ds, const UnitCell &ucell, const Parallel_Orbitals &pv, const ModuleSymmetry::Symmetry_rotation *p_symrot=nullptr)
Definition Exx_LRI.hpp:185
double post_process_Eexx(const double &Eexx_in) const
Definition Exx_LRI.hpp:248
std::pair< TA, TC > TAC
Definition Exx_LRI.h:57
int TA
Definition Exx_LRI.h:53
void cal_exx_ions(const UnitCell &ucell, const bool write_cv=false)
Definition Exx_LRI.hpp:70
void cal_exx_force(const int &nat)
Definition Exx_LRI.hpp:274
void cal_exx_stress(const double &omega, const double &lat0)
Definition Exx_LRI.hpp:297
Definition klist.h:13
Definition ORB_read.h:19
std::vector< double > cutoffs() const
Definition ORB_read.cpp:40
static void tick(const std::string &class_name_in, const std::string &name_in)
Use twice at a time: the first time, set start_flag to false; the second time, calculate the time dur...
Definition timer.cpp:57
Definition symmetry_rotation.h:16
std::map< int, std::map< std::pair< int, TC >, RI::Tensor< Tdata > > > restore_HR(const Symmetry &symm, const Atom *atoms, const Statistics &st, const char mode, const std::map< int, std::map< std::pair< int, TC >, RI::Tensor< Tdata > > > &HR_irreduceble) const
Definition symmetry_rotation_R.hpp:40
const std::map< Tap, std::set< TC > > & get_irreducible_sector() const
Definition symmetry_rotation.h:23
Definition parallel_orbitals.h:9
const Input_para & inp
Definition parameter.h:26
const System_para & globalv
Definition parameter.h:30
Definition unitcell.h:16
int *& iat2it
Definition unitcell.h:47
Atom * atoms
Definition unitcell.h:18
ModuleSymmetry::Symmetry symm
Definition unitcell.h:55
int & nat
Definition unitcell.h:46
ModuleBase::Vector3< double > & a2
Definition unitcell.h:36
ModuleBase::Vector3< double > & a3
Definition unitcell.h:36
int *& iat2ia
Definition unitcell.h:48
ModuleBase::Vector3< double > & a1
Definition unitcell.h:36
Statistics st
Definition unitcell.h:44
#define T
Definition exp.cpp:237
T cal_orbs_ccp(const T &orbs, const std::map< Conv_Coulomb_Pot_K::Coulomb_Type, std::vector< std::map< std::string, std::string > > > &coulomb_param, const double rmesh_times)
Exx_Info exx_info
Definition test_xc.cpp:29
int MY_RANK
global index of process
Definition global_variable.cpp:21
std::ofstream ofs_running
Definition global_variable.cpp:38
std::vector< std::array< T, N > > add(const std::vector< std::array< T, N > > &v1, const std::vector< std::array< T, N > > &v2)
Definition LRI_CV_Tools.hpp:219
std::map< int, std::map< int, std::map< Abfs::Vector3_Order< double >, RI::Tensor< Tdata > > > > get_CVws(const UnitCell &ucell, const std::map< TA, std::map< std::pair< TA, std::array< Tcell, 3 > >, RI::Tensor< Tdata > > > &CVs)
Definition LRI_CV_Tools.hpp:382
std::array< std::vector< T >, N > change_order(std::vector< std::array< T, N > > &&ds_in)
Definition LRI_CV_Tools.hpp:302
std::map< int, std::map< int, std::map< Abfs::Vector3_Order< double >, std::array< RI::Tensor< Tdata >, 3 > > > > get_dCVws(const UnitCell &ucell, const std::map< TA, std::map< std::pair< TA, std::array< Tcell, 3 > >, std::array< RI::Tensor< Tdata >, 3 > > > &dCVs)
Definition LRI_CV_Tools.hpp:409
void write_Cs_ao(const TLRI< T > &Vs, const std::string &file_path)
Definition write_ri_cv.hpp:50
std::array< std::array< std::map< TA, std::map< std::pair< TA, TC >, RI::Tensor< Tdata > > >, 3 >, 3 > cal_dMRs(const UnitCell &ucell, const std::array< std::map< TA, std::map< std::pair< TA, TC >, RI::Tensor< Tdata > > >, 3 > &dMs)
Definition LRI_CV_Tools.hpp:520
void write_Vs_abf(const TLRI< T > &Vs, const std::string &file_path)
Definition write_ri_cv.hpp:105
void TITLE(const std::string &class_name, const std::string &function_name, const bool disable)
Definition tool_title.cpp:18
std::vector< std::tuple< std::set< TA >, std::set< TA > > > get_2D_judge(const UnitCell &ucell, const Parallel_2D &pv)
Definition RI_2D_Comm.cpp:15
std::array< Tcell, 3 > Vector3_to_array3(const ModuleBase::Vector3< Tcell > &v)
Definition RI_Util.h:32
std::map< Conv_Coulomb_Pot_K::Coulomb_Method, std::pair< bool, std::map< Conv_Coulomb_Pot_K::Coulomb_Type, std::vector< std::map< std::string, std::string > > > > > update_coulomb_settings(const std::map< Conv_Coulomb_Pot_K::Coulomb_Type, std::vector< std::map< std::string, std::string > > > &coulomb_param, const UnitCell &ucell, const K_Vectors *p_kv)
Definition RI_Util.hpp:105
Parameter PARAM
Definition parameter.cpp:3
int abfs_Lmax
Definition exx_info.h:71
Exx_Info_RI info_ri
Definition exx_info.h:78
bool cal_force
calculate the force
Definition input_parameter.h:31
int nspin
LDA ; LSDA ; non-linear spin.
Definition input_parameter.h:84
bool cal_stress
calculate the stress
Definition input_parameter.h:32
std::string global_out_dir
global output directory
Definition system_parameter.h:42