ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
output_mulliken.h
Go to the documentation of this file.
1#ifndef OUTPUT_MULLIKEN_H
2#define OUTPUT_MULLIKEN_H
12
13#include <map>
14#include <vector>
15
16namespace ModuleIO
17{
18
20template <typename TK>
22{
23 public:
26 Output_DMK<TK>* output_dmk,
27 Parallel_Orbitals* ParaV,
28 CellIndex* cell_index,
29 const std::vector<int>& isk,
30 int nspin);
32 void write(int istep, std::string out_dir);
34 void print_atom_mag(const std::vector<std::vector<double>>& atom_chg, std::ostream& os);
36 std::vector<double> get_tot_chg();
38 std::vector<std::vector<double>> get_atom_chg();
40 std::map<std::vector<int>, double> get_orb_chg();
42 std::vector<std::vector<double>> get_atom_mulliken(std::vector<std::vector<double>>& atom_chg);
43
44 private:
45 /******************************************************************
46 * private functions
47 *******************************************************************/
49 void write_mulliken_nspin1(int istep,
50 const std::vector<double>& tot_chg,
51 const std::vector<std::vector<double>>& atom_chg,
52 std::map<std::vector<int>, double> orb_chg,
53 std::ofstream& os);
55 void write_mulliken_nspin2(int istep,
56 const std::vector<double>& tot_chg,
57 const std::vector<std::vector<double>>& atom_chg,
58 std::map<std::vector<int>, double> orb_chg,
59 std::ofstream& os);
61 void write_mulliken_nspin4(int istep,
62 const std::vector<double>& tot_chg,
63 const std::vector<std::vector<double>>& atom_chg,
64 std::map<std::vector<int>, double> orb_chg,
65 std::ofstream& os);
67 void set_nspin(int nspin_in);
69 void set_ParaV(Parallel_Orbitals* ParaV_in);
71 void collect_MW(ModuleBase::matrix& MecMulP, const ModuleBase::ComplexMatrix& mud, int nw, int isk);
74
75 private:
76 /******************************************************************
77 * private variables
78 *******************************************************************/
83 const std::vector<int>& isk_;
84 int nspin_;
86};
87
88template <typename TK>
90 hamilt::Hamilt<TK>* p_ham,
91 K_Vectors& kv,
93 const TwoCenterBundle& two_center_bundle,
94 const LCAO_Orbitals& orb,
95 UnitCell& ucell,
96 const Grid_Driver& gd,
97 const int istep,
98 const bool print)
99{
100 // 1) calculate and output Mulliken population charges and magnetic moments
101 if (PARAM.inp.out_mul)
102 {
103 auto cell_index
105 auto out_sk = ModuleIO::Output_Sk<TK>(p_ham, pv, PARAM.inp.nspin, kv.get_nks());
106 auto out_dmk = ModuleIO::Output_DMK<TK>(dynamic_cast<const elecstate::ElecStateLCAO<TK>*>(pelec)->get_DM(),
107 pv,
109 kv.get_nks());
110 auto mulp = ModuleIO::Output_Mulliken<TK>(&(out_sk), &(out_dmk), pv, &cell_index, kv.isk, PARAM.inp.nspin);
111 auto atom_chg = mulp.get_atom_chg();
113 ucell.atom_mulliken = mulp.get_atom_mulliken(atom_chg);
114 if (print && GlobalV::MY_RANK == 0)
115 {
117 cell_index.write_orb_info(PARAM.globalv.global_out_dir);
119 mulp.write(istep, PARAM.globalv.global_out_dir);
121 mulp.print_atom_mag(atom_chg, GlobalV::ofs_running);
122 }
123 }
124 // 2) calculate and output the magnetizations of each atom with projection method
125 if (PARAM.inp.onsite_radius > 0)
126 {
127 std::vector<std::vector<double>> atom_mag(ucell.nat, std::vector<double>(PARAM.inp.nspin, 0.0));
128 std::vector<ModuleBase::Vector3<int>> constrain(ucell.nat, ModuleBase::Vector3<int>(1, 1, 1));
130 = dynamic_cast<const elecstate::ElecStateLCAO<TK>*>(pelec)->get_DM()->get_DMR_pointer(1);
131 std::vector<double> moments;
132 std::vector<double> mag_x(ucell.nat, 0.0);
133 std::vector<double> mag_y(ucell.nat, 0.0);
134 std::vector<double> mag_z(ucell.nat, 0.0);
135 auto atomLabels = ucell.get_atomLabels();
136 if(PARAM.inp.nspin == 2)
137 {
138 auto sc_lambda
140 kv.kvec_d,
141 dynamic_cast<hamilt::HamiltLCAO<TK, double>*>(p_ham)->getHR(),
142 ucell,
143 &gd,
144 two_center_bundle.overlap_orb_onsite.get(),
145 orb.cutoffs());
146 dynamic_cast<const elecstate::ElecStateLCAO<TK>*>(pelec)->get_DM()->switch_dmr(2);
147 moments = sc_lambda->cal_moment(dmr, constrain);
148 dynamic_cast<const elecstate::ElecStateLCAO<TK>*>(pelec)->get_DM()->switch_dmr(0);
149 delete sc_lambda;
150 //const std::vector<std::string> title = {"Total Magnetism (uB)", ""};
151 //const std::vector<std::string> fmts = {"%-26s", "%20.10f"};
152 //FmtTable table(title, ucell.nat, fmts, {FmtTable::Align::RIGHT, FmtTable::Align::LEFT});
153 for(int iat=0;iat<ucell.nat;iat++)
154 {
155 atom_mag[iat][0] = 0.0;
156 atom_mag[iat][1] = moments[iat];
157 // mag_z[iat] = moments[iat];
158 }
159 //table << atomLabels << mag_z;
160 //GlobalV::ofs_running << table.str() << std::endl;
161 }
162 else if(PARAM.inp.nspin == 4)
163 {
164 auto sc_lambda = new hamilt::DeltaSpin<hamilt::OperatorLCAO<std::complex<double>, std::complex<double>>>(
165 nullptr,
166 kv.kvec_d,
167 dynamic_cast<hamilt::HamiltLCAO<std::complex<double>, std::complex<double>>*>(p_ham)->getHR(),
168 ucell,
169 &gd,
170 two_center_bundle.overlap_orb_onsite.get(),
171 orb.cutoffs());
172 moments = sc_lambda->cal_moment(dmr, constrain);
173 delete sc_lambda;
174 //const std::vector<std::string> title = {"Total Magnetism (uB)", "", "", ""};
175 //const std::vector<std::string> fmts = {"%-26s", "%20.10f", "%20.10f", "%20.10f"};
176 //FmtTable table(title, ucell.nat, fmts, {FmtTable::Align::RIGHT, FmtTable::Align::LEFT});
177 for(int iat=0;iat<ucell.nat;iat++)
178 {
179 atom_mag[iat][0] = 0.0;
180 atom_mag[iat][1] = moments[iat*3];
181 atom_mag[iat][2] = moments[iat*3+1];
182 atom_mag[iat][3] = moments[iat*3+2];
183 //mag_x[iat] = moments[iat*3];
184 //mag_y[iat] = moments[iat*3+1];
185 //mag_z[iat] = moments[iat*3+2];
186 }
187 //table << atomLabels << mag_x << mag_y << mag_z;
188 //GlobalV::ofs_running << table.str() << std::endl;
189 }
190 ucell.atom_mulliken = atom_mag;
191 }
192}
193
194} // namespace ModuleIO
195
196#endif // OUTPUT_MULLIKEN_H
This class is used to get the information of the atoms and orbitals indices in the unit cell.
Definition cell_index.h:21
Definition sltk_grid_driver.h:43
Definition klist.h:13
int get_nks() const
Definition klist.h:68
std::vector< ModuleBase::Vector3< double > > kvec_d
Cartesian coordinates of k points.
Definition klist.h:16
std::vector< int > isk
ngk, number of plane waves for each k point
Definition klist.h:21
Definition ORB_read.h:19
std::vector< double > cutoffs() const
Definition ORB_read.cpp:40
Definition complexmatrix.h:14
3 elements vector
Definition vector3.h:22
Definition matrix.h:19
Definition output_dmk.h:11
the output interface to write the Mulliken population charges
Definition output_mulliken.h:22
CellIndex * cell_index_
Definition output_mulliken.h:82
std::vector< double > get_tot_chg()
get total charge
Definition output_mulliken.cpp:305
void print_atom_mag(const std::vector< std::vector< double > > &atom_chg, std::ostream &os)
print atom mag to running log file
Definition output_mulliken.cpp:432
void set_ParaV(Parallel_Orbitals *ParaV_in)
set orbital parallel info
Definition output_mulliken.cpp:294
std::vector< std::vector< double > > get_atom_chg()
get atom charge
Definition output_mulliken.cpp:321
int nspin_
Definition output_mulliken.h:84
const std::vector< int > & isk_
Definition output_mulliken.h:83
Output_DMK< TK > * output_dmk_
Definition output_mulliken.h:80
void write_mulliken_nspin2(int istep, const std::vector< double > &tot_chg, const std::vector< std::vector< double > > &atom_chg, std::map< std::vector< int >, double > orb_chg, std::ofstream &os)
write mulliken.txt for the case of nspin=2
Definition output_mulliken.cpp:128
std::map< std::vector< int >, double > get_orb_chg()
get orbital charge
Definition output_mulliken.cpp:342
Output_Sk< TK > * output_sk_
Definition output_mulliken.h:79
void write_mulliken_nspin1(int istep, const std::vector< double > &tot_chg, const std::vector< std::vector< double > > &atom_chg, std::map< std::vector< int >, double > orb_chg, std::ofstream &os)
write mulliken.txt for the case of nspin=1
Definition output_mulliken.cpp:65
void collect_MW(ModuleBase::matrix &MecMulP, const ModuleBase::ComplexMatrix &mud, int nw, int isk)
collect_mw from matrix multiplication result
Definition output_mulliken.cpp:371
void write_mulliken_nspin4(int istep, const std::vector< double > &tot_chg, const std::vector< std::vector< double > > &atom_chg, std::map< std::vector< int >, double > orb_chg, std::ofstream &os)
write mulliken.txt for the case of nspin=4
Definition output_mulliken.cpp:208
Parallel_Orbitals * ParaV_
Definition output_mulliken.h:81
std::vector< std::vector< double > > get_atom_mulliken(std::vector< std::vector< double > > &atom_chg)
returun atom_mulliken for updateing STRU file
Definition output_mulliken.cpp:502
void set_nspin(int nspin_in)
set nspin
Definition output_mulliken.cpp:283
void write(int istep, std::string out_dir)
the outer interface to write the Mulliken population charges
Definition output_mulliken.cpp:29
void cal_orbMulP()
mulliken population = trace(dm*overlap)
ModuleBase::matrix orbMulP_
Definition output_mulliken.h:85
Definition output_sk.h:12
Definition parallel_orbitals.h:9
const Input_para & inp
Definition parameter.h:26
const System_para & globalv
Definition parameter.h:30
Definition two_center_bundle.h:11
std::unique_ptr< TwoCenterIntegrator > overlap_orb_onsite
Definition two_center_bundle.h:46
Definition unitcell.h:16
std::vector< int > get_atomCounts() const
get atomCounts, which is a vector of element type with atom number
Definition unitcell.cpp:120
std::vector< std::string > get_atomLabels() const
get atom labels
Definition unitcell.cpp:112
std::vector< std::vector< double > > atom_mulliken
Definition unitcell.h:22
int & nat
Definition unitcell.h:46
std::vector< std::vector< int > > get_lnchiCounts() const
get lnchiCounts, which is a vector of element type with the l:nchi vector
Definition unitcell.cpp:128
Definition elecstate_lcao.h:15
Definition elecstate.h:15
Definition dspin_lcao.h:20
Definition hcontainer.h:144
Definition hamilt_lcao.h:33
Definition hamilt.h:16
a new formatter library for formatting data
int MY_RANK
global index of process
Definition global_variable.cpp:21
std::ofstream ofs_running
Definition global_variable.cpp:38
This class has two functions: restart psi from the previous calculation, and write psi to the disk.
Definition cal_dos.h:9
void cal_mag(Parallel_Orbitals *pv, hamilt::Hamilt< TK > *p_ham, K_Vectors &kv, elecstate::ElecState *pelec, const TwoCenterBundle &two_center_bundle, const LCAO_Orbitals &orb, UnitCell &ucell, const Grid_Driver &gd, const int istep, const bool print)
Definition output_mulliken.h:89
Parameter PARAM
Definition parameter.cpp:3
bool out_mul
qifeng add 2019-9-10
Definition input_parameter.h:373
double onsite_radius
radius of the sphere for onsite projection (Bohr)
Definition input_parameter.h:560
int nspin
LDA ; LSDA ; non-linear spin.
Definition input_parameter.h:84
std::string global_out_dir
global output directory
Definition system_parameter.h:42