ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
get_pchg_pw.h
Go to the documentation of this file.
1#ifndef GET_PCHG_PW_H
2#define GET_PCHG_PW_H
3
4#include "source_io/cube_io.h"
6
7namespace ModuleIO
8{
9template <typename Device>
10void get_pchg_pw(const std::vector<int>& out_pchg,
11 const int nbands,
12 const int nspin,
13 const int nxyz,
14 const int chr_ngmc,
15 UnitCell* ucell,
16 const psi::Psi<std::complex<double>, Device>* kspw_psi,
17 const ModulePW::PW_Basis* pw_rhod,
18 const ModulePW::PW_Basis_K* pw_wfc,
19 const Device* ctx,
20 const Parallel_Grid& pgrid,
21 const std::string& global_out_dir,
22 const bool if_separate_k,
23 const K_Vectors& kv,
24 const int kpar,
25 const int my_pool,
26 const Charge* chr) // Charge class is needed for the charge density reduce
27{
28 // Get necessary parameters from kv
29 const int nks = kv.get_nks(); // current process pool k-point count
30 const int nkstot = kv.get_nkstot(); // total k-point count
31
32 // Loop over k-parallelism
33 for (int ip = 0; ip < kpar; ++ip)
34 {
35 if (my_pool != ip)
36 {
37 continue;
38 }
39
40 // bands_picked is a vector of 0s and 1s, where 1 means the band is picked to output
41 std::vector<int> bands_picked(nbands, 0);
42
43 // Check if length of out_pchg is valid
44 if (static_cast<int>(out_pchg.size()) > nbands)
45 {
46 ModuleBase::WARNING_QUIT("ModuleIO::get_pchg_pw",
47 "The number of bands specified by `out_pchg` in the "
48 "INPUT file exceeds `nbands`!");
49 }
50
51 // Check if all elements in bands_picked are 0 or 1
52 for (int value: out_pchg)
53 {
54 if (value != 0 && value != 1)
55 {
56 ModuleBase::WARNING_QUIT("ModuleIO::get_pchg_pw",
57 "The elements of `out_pchg` must be either 0 or 1. "
58 "Invalid values found!");
59 }
60 }
61
62 // Fill bands_picked with values from out_pchg
63 // Remaining bands are already set to 0
64 int length = std::min(static_cast<int>(out_pchg.size()), nbands);
65 for (int i = 0; i < length; ++i)
66 {
67 // out_pchg rely on function parse_expression
68 bands_picked[i] = static_cast<int>(out_pchg[i]);
69 }
70
71 // Allocate host memory
72 std::vector<std::complex<double>> wfcr(nxyz);
73 std::vector<std::vector<double>> rho_band(nspin, std::vector<double>(nxyz));
74
75 // Allocate device memory
76 std::complex<double>* wfcr_device = nullptr;
77 if (!std::is_same<Device, base_device::DEVICE_CPU>::value)
78 {
80 }
81
82 for (int ib = 0; ib < nbands; ++ib)
83 {
84 // Skip the loop iteration if bands_picked[ib] is 0
85 if (!bands_picked[ib])
86 {
87 continue;
88 }
89
90 for (int is = 0; is < nspin; ++is)
91 {
92 std::fill(rho_band[is].begin(), rho_band[is].end(), 0.0);
93 }
94
95 if (if_separate_k)
96 {
97 for (int ik = 0; ik < nks; ++ik)
98 {
99 const int ikstot = kv.ik2iktot[ik]; // global k-point index
100 const int spin_index = kv.isk[ik]; // spin index
101 const int k_number = ikstot % (nkstot / nspin) + 1; // k-point number, starting from 1
102
103 kspw_psi->fix_k(ik);
104
105 // FFT on device and copy result back to host
106 if (std::is_same<Device, base_device::DEVICE_CPU>::value)
107 {
108 pw_wfc->recip_to_real(ctx, &kspw_psi[0](ib, 0), wfcr.data(), ik);
109 }
110 else
111 {
112 pw_wfc->recip_to_real(ctx, &kspw_psi[0](ib, 0), wfcr_device, ik);
113
115 base_device::DEVICE_CPU,
116 Device>()(wfcr.data(), wfcr_device, nxyz);
117 }
118
119 // To ensure the normalization of charge density in multi-k calculation (if if_separate_k is true)
120 double wg_sum_k = 0.0;
121 if (nspin == 1)
122 {
123 wg_sum_k = 2.0;
124 }
125 else if (nspin == 2)
126 {
127 wg_sum_k = 1.0;
128 }
129 else
130 {
131 ModuleBase::WARNING_QUIT("ModuleIO::get_pchg_pw",
132 "Real space partial charge output currently do not support "
133 "noncollinear polarized calculation (nspin = 4)!");
134 }
135
136 double w1 = static_cast<double>(wg_sum_k / ucell->omega);
137
138 for (int i = 0; i < nxyz; ++i)
139 {
140 rho_band[spin_index][i] = std::norm(wfcr[i]) * w1;
141 }
142
143 std::stringstream ssc;
144 ssc << global_out_dir << "pchgi" << ib + 1 << "s" << spin_index + 1 << "k" << k_number << ".cube";
145
147 rho_band[spin_index].data(),
148 spin_index,
149 nspin,
150 0,
151 ssc.str(),
152 0.0,
153 ucell,
154 11,
155 1,
156 true); // reduce_all_pool is true
157 }
158 }
159 else
160 {
161 for (int ik = 0; ik < nks; ++ik)
162 {
163 const int ikstot = kv.ik2iktot[ik]; // global k-point index
164 const int spin_index = kv.isk[ik]; // spin index
165 const int k_number = ikstot % (nkstot / nspin) + 1; // k-point number, starting from 1
166
167 kspw_psi->fix_k(ik);
168
169 // FFT on device and copy result back to host
170 if (std::is_same<Device, base_device::DEVICE_CPU>::value)
171 {
172 pw_wfc->recip_to_real(ctx, &kspw_psi[0](ib, 0), wfcr.data(), ik);
173 }
174 else
175 {
176 pw_wfc->recip_to_real(ctx, &kspw_psi[0](ib, 0), wfcr_device, ik);
177
179 base_device::DEVICE_CPU,
180 Device>()(wfcr.data(), wfcr_device, nxyz);
181 }
182
183 double w1 = static_cast<double>(kv.wk[ik] / ucell->omega);
184
185 for (int i = 0; i < nxyz; ++i)
186 {
187 rho_band[spin_index][i] += std::norm(wfcr[i]) * w1;
188 }
189 }
190
191#ifdef __MPI
192 // Reduce the charge density across all pools if kpar > 1
193 if (kpar > 1 && chr != nullptr)
194 {
195 for (int is = 0; is < nspin; ++is)
196 {
197 chr->reduce_diff_pools(rho_band[is].data());
198 }
199 }
200#endif
201
202 // Symmetrize the charge density, otherwise the results are incorrect if the symmetry is on
203 // std::cout << " Symmetrizing band-decomposed charge density..." << std::endl;
204 Symmetry_rho srho;
205 for (int is = 0; is < nspin; ++is)
206 {
207 // Use vector instead of raw pointers
208 std::vector<double*> rho_save_pointers(nspin);
209 for (int s = 0; s < nspin; ++s)
210 {
211 rho_save_pointers[s] = rho_band[s].data();
212 }
213
214 std::vector<std::vector<std::complex<double>>> rhog(nspin,
215 std::vector<std::complex<double>>(chr_ngmc));
216
217 // Convert vector of vectors to vector of pointers
218 std::vector<std::complex<double>*> rhog_pointers(nspin);
219 for (int s = 0; s < nspin; ++s)
220 {
221 rhog_pointers[s] = rhog[s].data();
222 }
223
224 srho.begin(is,
225 rho_save_pointers.data(),
226 rhog_pointers.data(),
227 chr_ngmc,
228 nullptr,
229 pw_rhod,
230 ucell->symm);
231 }
232
233 for (int is = 0; is < nspin; ++is)
234 {
235 std::stringstream ssc;
236 ssc << global_out_dir << "pchgi" << ib + 1 << "s" << is + 1 << ".cube";
237
238 ModuleIO::write_vdata_palgrid(pgrid, rho_band[is].data(), is, nspin, 0, ssc.str(), 0.0, ucell);
239 }
240 } // else if_separate_k is false
241 } // end of ib loop over nbands
242 } // end of ip loop over kpar
243} // get_pchg_pw
244} // namespace ModuleIO
245
246#endif // GET_PCHG_PW_H
Definition charge.h:18
void reduce_diff_pools(double *array_rho) const
Reduce among different pools If NPROC_IN_POOLs are all the same, use GlobalV::KP_WORLD else,...
Definition charge_mpi.cpp:27
Definition klist.h:13
int get_nks() const
Definition klist.h:68
std::vector< int > ik2iktot
[nks] map ik to the global index of k points
Definition klist.h:128
std::vector< int > isk
ngk, number of plane waves for each k point
Definition klist.h:21
int get_nkstot() const
Definition klist.h:73
std::vector< double > wk
Direct coordinates of k points.
Definition klist.h:18
Special pw_basis class. It includes different k-points.
Definition pw_basis_k.h:57
void recip_to_real(const Device *ctx, const std::complex< FPTYPE > *in, std::complex< FPTYPE > *out, const int ik, const bool add=false, const FPTYPE factor=1.0) const
Definition xc3_mock.h:94
A class which can convert a function of "r" to the corresponding linear superposition of plane waves ...
Definition pw_basis.h:56
Definition parallel_grid.h:8
Definition symmetry_rho.h:9
void begin(const int &spin_now, const Charge &CHR, const ModulePW::PW_Basis *pw, ModuleSymmetry::Symmetry &symm) const
Definition symmetry_rho.cpp:13
Definition unitcell.h:17
ModuleSymmetry::Symmetry symm
Definition unitcell.h:57
double & omega
Definition unitcell.h:34
Definition psi.h:37
void WARNING_QUIT(const std::string &, const std::string &)
Combine the functions of WARNING and QUIT.
Definition test_delley.cpp:14
Definition cal_dos.h:9
void get_pchg_pw(const std::vector< int > &out_pchg, const int nbands, const int nspin, const int nxyz, const int chr_ngmc, UnitCell *ucell, const psi::Psi< std::complex< double >, Device > *kspw_psi, const ModulePW::PW_Basis *pw_rhod, const ModulePW::PW_Basis_K *pw_wfc, const Device *ctx, const Parallel_Grid &pgrid, const std::string &global_out_dir, const bool if_separate_k, const K_Vectors &kv, const int kpar, const int my_pool, const Charge *chr)
Definition get_pchg_pw.h:10
void write_vdata_palgrid(const Parallel_Grid &pgrid, const double *const data, const int is, const int nspin, const int iter, const std::string &fn, const double ef, const UnitCell *const ucell, const int precision=11, const int out_fermi=1, const bool reduce_all_pool=false)
write volumetric data on the parallized grid into a .cube file
Definition write_cube.cpp:15
Definition memory_op.h:17