ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
spar_hsr.h
Go to the documentation of this file.
1#ifndef SPARSE_FORMAT_HSR_H
2#define SPARSE_FORMAT_HSR_H
3
6#include "source_lcao/module_dftu/dftu.h" // mohan add 20251107
7
8namespace sparse_format
9{
10#ifdef __MPI
11// Synchronize all processes' R coordinates,
12// otherwise HS_Arrays.all_R_coor would have different sizes on different processes,
13// causing MPI communication errors.
14void sync_all_R_coor(std::set<Abfs::Vector3_Order<int>>& all_R_coor, MPI_Comm comm);
15#endif
16
17template <typename T>
18std::set<Abfs::Vector3_Order<int>> get_R_range(const hamilt::HContainer<T>& hR)
19{
20 std::set<Abfs::Vector3_Order<int>> all_R_coor;
21 for (int iap = 0; iap < hR.size_atom_pairs(); ++iap)
22 {
23 const hamilt::AtomPair<T>& atom_pair = hR.get_atom_pair(iap);
24 for (int iR = 0; iR < atom_pair.get_R_size(); ++iR)
25 {
26 const auto& r_index = atom_pair.get_R_index(iR);
27 Abfs::Vector3_Order<int> dR(r_index.x, r_index.y, r_index.z);
28 all_R_coor.insert(dR);
29 }
30 }
31
32#ifdef __MPI
33 // Fix: Sync all_R_coor across processes
34 sparse_format::sync_all_R_coor(all_R_coor, MPI_COMM_WORLD);
35#endif
36
37 return all_R_coor;
38};
39
40using TAC = std::pair<int, std::array<int, 3>>;
41template <typename TK>
42void cal_HSR(const UnitCell& ucell,
43 Plus_U &dftu, // mohan add 2025-11-07
44 const Parallel_Orbitals& pv,
45 LCAO_HS_Arrays& HS_Arrays,
46 const Grid_Driver& grid,
47 const int& current_spin,
48 const double& sparse_thr,
49 const int (&nmp)[3],
51#ifdef __EXX
52 ,
53 const std::vector<std::map<int, std::map<TAC, RI::Tensor<double>>>>* Hexxd = nullptr,
54 const std::vector<std::map<int, std::map<TAC, RI::Tensor<std::complex<double>>>>>* Hexxc = nullptr
55#endif
56 );
57
58template <typename TI, typename TO = TI>
60 const double& sparse_thr,
61 const hamilt::HContainer<TI>& hR,
62 std::map<Abfs::Vector3_Order<int>, std::map<size_t, std::map<size_t, TO>>>& target);
63
64void clear_zero_elements(LCAO_HS_Arrays& HS_Arrays, const int& current_spin, const double& sparse_thr);
65
66void destroy_HS_R_sparse(LCAO_HS_Arrays& HS_Arrays);
67
68} // namespace sparse_format
69
70#endif
Definition abfs-vector3_order.h:16
Definition sltk_grid_driver.h:43
Definition LCAO_HS_arrays.hpp:9
Definition parallel_orbitals.h:9
Definition dftu.h:20
Definition unitcell.h:17
Definition atom_pair.h:42
ModuleBase::Vector3< int > get_R_index(const int &index) const
Definition atom_pair.cpp:795
size_t get_R_size() const
Definition atom_pair.h:288
Definition hcontainer.h:144
size_t size_atom_pairs() const
calculate number of AtomPairs for current R index
Definition hcontainer.cpp:576
AtomPair< T > & get_atom_pair(int i, int j) const
return a reference of AtomPair with index of atom I and J in atom_pairs
Definition hcontainer.cpp:353
Definition hamilt.h:16
Definition spar_dh.h:13
std::pair< int, std::array< int, 3 > > TAC
Definition spar_hsr.h:40
void destroy_HS_R_sparse(LCAO_HS_Arrays &HS_Arrays)
Definition spar_hsr.cpp:345
void cal_HSR(const UnitCell &ucell, Plus_U &dftu, const Parallel_Orbitals &pv, LCAO_HS_Arrays &HS_Arrays, const Grid_Driver &grid, const int &current_spin, const double &sparse_thr, const int(&nmp)[3], hamilt::Hamilt< TK > *p_ham)
Definition spar_hsr.cpp:71
std::set< Abfs::Vector3_Order< int > > get_R_range(const hamilt::HContainer< T > &hR)
Definition spar_hsr.h:18
void clear_zero_elements(LCAO_HS_Arrays &HS_Arrays, const int &current_spin, const double &sparse_thr)
Definition spar_hsr.cpp:231
void sync_all_R_coor(std::set< Abfs::Vector3_Order< int > > &all_R_coor, MPI_Comm comm)
Definition spar_hsr.cpp:11
void cal_HContainer(const Parallel_Orbitals &pv, const double &sparse_thr, const hamilt::HContainer< TI > &hR, std::map< Abfs::Vector3_Order< int >, std::map< size_t, std::map< size_t, TO > > > &target)
Definition spar_hsr.cpp:188
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