Skip to content

Commit ec23da3

Browse files
committed
use original format
1 parent e2bb480 commit ec23da3

File tree

3 files changed

+47
-61
lines changed

3 files changed

+47
-61
lines changed

source/module_esolver/esolver_ks_lcao.cpp

Lines changed: 1 addition & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -1153,13 +1153,7 @@ void ESolver_KS_LCAO<TK, TR>::after_scf(const int istep)
11531153
// mulliken charge analysis
11541154
if (PARAM.inp.out_mul)
11551155
{
1156-
auto mulp = ModuleIO::Output_Mulliken<TK>(&(this->pv),
1157-
this->p_hamilt,
1158-
this->kv,
1159-
dynamic_cast<const elecstate::ElecStateLCAO<TK>*>(this->pelec),
1160-
GlobalC::ucell,
1161-
PARAM.inp.nspin);
1162-
mulp.cal_mag(GlobalC::ucell, istep, true);
1156+
ModuleIO::cal_mag(&(this->pv), this->p_hamilt, this->kv, this->pelec, GlobalC::ucell, istep, true);
11631157
}
11641158
}
11651159

source/module_io/output_mulliken.cpp

Lines changed: 8 additions & 46 deletions
Original file line numberDiff line numberDiff line change
@@ -12,20 +12,16 @@ namespace ModuleIO
1212
{
1313

1414
template <typename TK>
15-
Output_Mulliken<TK>::Output_Mulliken(Parallel_Orbitals* pv,
16-
hamilt::Hamilt<TK>* p_ham,
17-
const K_Vectors& kv,
18-
const elecstate::ElecStateLCAO<TK>* pelec,
19-
const UnitCell& ucell,
20-
const int nspin)
15+
Output_Mulliken<TK>::Output_Mulliken(Output_Sk<TK>* output_sk,
16+
Output_DMK<TK>* output_dmk,
17+
Parallel_Orbitals* ParaV,
18+
CellIndex* cell_index,
19+
const std::vector<int>& isk,
20+
int nspin)
21+
: output_sk_(output_sk), output_dmk_(output_dmk), ParaV_(ParaV), cell_index_(cell_index), isk_(isk), nspin_(nspin)
2122
{
22-
this->cell_index_
23-
= new CellIndex(ucell.get_atomLabels(), ucell.get_atomCounts(), ucell.get_lnchiCounts(), PARAM.inp.nspin);
24-
this->output_sk_ = new ModuleIO::Output_Sk<TK>(p_ham, pv, PARAM.inp.nspin, kv.get_nks());
25-
this->output_dmk_ = new ModuleIO::Output_DMK<TK>(pelec->get_DM(), pv, PARAM.inp.nspin, kv.get_nks());
26-
this->isk_ = kv.isk;
2723
this->set_nspin(nspin);
28-
this->set_ParaV(pv);
24+
this->set_ParaV(ParaV);
2925
this->cal_orbMulP();
3026
}
3127

@@ -65,40 +61,6 @@ void Output_Mulliken<TK>::write(int istep, std::string out_dir)
6561
os.close();
6662
}
6763

68-
template <>
69-
void Output_Mulliken<std::complex<double>>::cal_mag(UnitCell& ucell, const int istep, const bool print)
70-
{
71-
auto atom_chg = this->get_atom_chg();
72-
/// used in updating mag info in STRU file
73-
ucell.atom_mulliken = this->get_atom_mulliken(atom_chg);
74-
if (print && GlobalV::MY_RANK == 0)
75-
{
76-
/// write the Orbital file
77-
this->cell_index_->write_orb_info(PARAM.globalv.global_out_dir);
78-
/// write mulliken.txt
79-
this->write(istep, PARAM.globalv.global_out_dir);
80-
/// write atomic mag info in running log file
81-
this->print_atom_mag(atom_chg, GlobalV::ofs_running);
82-
}
83-
}
84-
85-
template <>
86-
void Output_Mulliken<double>::cal_mag(UnitCell& ucell, const int istep, const bool print)
87-
{
88-
auto atom_chg = this->get_atom_chg();
89-
/// used in updating mag info in STRU file
90-
ucell.atom_mulliken = this->get_atom_mulliken(atom_chg);
91-
if (print && GlobalV::MY_RANK == 0)
92-
{
93-
/// write the Orbital file
94-
this->cell_index_->write_orb_info(PARAM.globalv.global_out_dir);
95-
/// write mulliken.txt
96-
this->write(istep, PARAM.globalv.global_out_dir);
97-
/// write atomic mag info in running log file
98-
this->print_atom_mag(atom_chg, GlobalV::ofs_running);
99-
}
100-
}
101-
10264
template <typename TK>
10365
void Output_Mulliken<TK>::write_mulliken_nspin1(int istep,
10466
const std::vector<double>& tot_chg,

source/module_io/output_mulliken.h

Lines changed: 38 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -20,13 +20,12 @@ class Output_Mulliken
2020
{
2121
public:
2222
/// constructor of Output_Mulliken
23-
Output_Mulliken(Parallel_Orbitals* pv,
24-
hamilt::Hamilt<TK>* p_ham,
25-
const K_Vectors& kv,
26-
const elecstate::ElecStateLCAO<TK>* pelec,
27-
const UnitCell& ucell,
28-
const int nspin);
29-
void cal_mag(UnitCell& ucell, const int istep, const bool print);
23+
Output_Mulliken(Output_Sk<TK>* output_sk,
24+
Output_DMK<TK>* output_dmk,
25+
Parallel_Orbitals* ParaV,
26+
CellIndex* cell_index,
27+
const std::vector<int>& isk,
28+
int nspin);
3029
/// the outer interface to write the Mulliken population charges
3130
void write(int istep, std::string out_dir);
3231
/// print atom mag to running log file
@@ -79,11 +78,42 @@ class Output_Mulliken
7978
Output_DMK<TK>* output_dmk_ = nullptr;
8079
Parallel_Orbitals* ParaV_ = nullptr;
8180
CellIndex* cell_index_ = nullptr;
82-
std::vector<int> isk_;
81+
const std::vector<int>& isk_;
8382
int nspin_;
8483
ModuleBase::matrix orbMulP_;
8584
};
8685

86+
template <typename TK>
87+
void cal_mag(Parallel_Orbitals* pv,
88+
hamilt::Hamilt<TK>* p_ham,
89+
K_Vectors& kv,
90+
elecstate::ElecState* pelec,
91+
UnitCell& ucell,
92+
const int istep,
93+
const bool print)
94+
{
95+
auto cell_index
96+
= CellIndex(ucell.get_atomLabels(), ucell.get_atomCounts(), ucell.get_lnchiCounts(), PARAM.inp.nspin);
97+
auto out_sk = ModuleIO::Output_Sk<TK>(p_ham, pv, PARAM.inp.nspin, kv.get_nks());
98+
auto out_dmk = ModuleIO::Output_DMK<TK>(dynamic_cast<const elecstate::ElecStateLCAO<TK>*>(pelec)->get_DM(),
99+
pv,
100+
PARAM.inp.nspin,
101+
kv.get_nks());
102+
auto mulp = ModuleIO::Output_Mulliken<TK>(&(out_sk), &(out_dmk), pv, &cell_index, kv.isk, PARAM.inp.nspin);
103+
auto atom_chg = mulp.get_atom_chg();
104+
/// used in updating mag info in STRU file
105+
ucell.atom_mulliken = mulp.get_atom_mulliken(atom_chg);
106+
if (print && GlobalV::MY_RANK == 0)
107+
{
108+
/// write the Orbital file
109+
cell_index.write_orb_info(PARAM.globalv.global_out_dir);
110+
/// write mulliken.txt
111+
mulp.write(istep, PARAM.globalv.global_out_dir);
112+
/// write atomic mag info in running log file
113+
mulp.print_atom_mag(atom_chg, GlobalV::ofs_running);
114+
}
115+
}
116+
87117
} // namespace ModuleIO
88118

89119
#endif // OUTPUT_MULLIKEN_H

0 commit comments

Comments
 (0)