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
5
6namespace sparse_format
7{
8#ifdef __MPI
9// Synchronize all processes' R coordinates,
10// otherwise HS_Arrays.all_R_coor would have different sizes on different processes,
11// causing MPI communication errors.
12void sync_all_R_coor(std::set<Abfs::Vector3_Order<int>>& all_R_coor, MPI_Comm comm);
13#endif
14
15template <typename T>
16std::set<Abfs::Vector3_Order<int>> get_R_range(const hamilt::HContainer<T>& hR)
17{
18 std::set<Abfs::Vector3_Order<int>> all_R_coor;
19 for (int iap = 0; iap < hR.size_atom_pairs(); ++iap)
20 {
21 const hamilt::AtomPair<T>& atom_pair = hR.get_atom_pair(iap);
22 for (int iR = 0; iR < atom_pair.get_R_size(); ++iR)
23 {
24 const auto& r_index = atom_pair.get_R_index(iR);
25 Abfs::Vector3_Order<int> dR(r_index.x, r_index.y, r_index.z);
26 all_R_coor.insert(dR);
27 }
28 }
29
30#ifdef __MPI
31 // Fix: Sync all_R_coor across processes
32 sparse_format::sync_all_R_coor(all_R_coor, MPI_COMM_WORLD);
33#endif
34
35 return all_R_coor;
36};
37
38using TAC = std::pair<int, std::array<int, 3>>;
39void cal_HSR(const UnitCell& ucell,
40 const Parallel_Orbitals& pv,
41 LCAO_HS_Arrays& HS_Arrays,
42 const Grid_Driver& grid,
43 const int& current_spin,
44 const double& sparse_thr,
45 const int (&nmp)[3],
46 hamilt::Hamilt<std::complex<double>>* p_ham
47#ifdef __EXX
48 ,
49 const std::vector<std::map<int, std::map<TAC, RI::Tensor<double>>>>* Hexxd = nullptr,
50 const std::vector<std::map<int, std::map<TAC, RI::Tensor<std::complex<double>>>>>* Hexxc = nullptr
51#endif
52);
53
55 const int& current_spin,
56 const double& sparse_threshold,
58 std::map<Abfs::Vector3_Order<int>, std::map<size_t, std::map<size_t, double>>>& target);
59
61 const Parallel_Orbitals& pv,
62 const int& current_spin,
63 const double& sparse_threshold,
64 const hamilt::HContainer<std::complex<double>>& hR,
65 std::map<Abfs::Vector3_Order<int>, std::map<size_t, std::map<size_t, std::complex<double>>>>& target);
66
68 const Parallel_Orbitals& pv,
69 const int& current_spin,
70 const double& sparse_threshold,
72 std::map<Abfs::Vector3_Order<int>, std::map<size_t, std::map<size_t, std::complex<double>>>>& target);
73
74void clear_zero_elements(LCAO_HS_Arrays& HS_Arrays, const int& current_spin, const double& sparse_thr);
75
76void destroy_HS_R_sparse(LCAO_HS_Arrays& HS_Arrays);
77
78} // namespace sparse_format
79
80#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 unitcell.h:16
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:38
void destroy_HS_R_sparse(LCAO_HS_Arrays &HS_Arrays)
Definition spar_hsr.cpp:426
void cal_HSR(const UnitCell &ucell, 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< std::complex< double > > *p_ham)
Definition spar_hsr.cpp:70
std::set< Abfs::Vector3_Order< int > > get_R_range(const hamilt::HContainer< T > &hR)
Definition spar_hsr.h:16
void cal_HContainer_d(const Parallel_Orbitals &pv, const int &current_spin, const double &sparse_threshold, const hamilt::HContainer< double > &hR, std::map< Abfs::Vector3_Order< int >, std::map< size_t, std::map< size_t, double > > > &target)
Definition spar_hsr.cpp:182
void cal_HContainer_td(const Parallel_Orbitals &pv, const int &current_spin, const double &sparse_threshold, const hamilt::HContainer< double > &hR, std::map< Abfs::Vector3_Order< int >, std::map< size_t, std::map< size_t, std::complex< double > > > > &target)
Definition spar_hsr.cpp:268
void clear_zero_elements(LCAO_HS_Arrays &HS_Arrays, const int &current_spin, const double &sparse_thr)
Definition spar_hsr.cpp:312
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_cd(const Parallel_Orbitals &pv, const int &current_spin, const double &sparse_threshold, const hamilt::HContainer< std::complex< double > > &hR, std::map< Abfs::Vector3_Order< int >, std::map< size_t, std::map< size_t, std::complex< double > > > > &target)
Definition spar_hsr.cpp:225
base device SOURCES math_dngvd_test cpp endif() if(ENABLE_GOOGLEBENCH) AddTest(TARGET PERF_MODULE_HSOLVER_KERNELS LIBS parameter $
Definition CMakeLists.txt:10