Skip to content

Commit

Permalink
refactor: update kpoints format (#1474)
Browse files Browse the repository at this point in the history
  • Loading branch information
hongriTianqi committed Nov 4, 2022
1 parent 84324f4 commit f145042
Show file tree
Hide file tree
Showing 2 changed files with 30 additions and 20 deletions.
46 changes: 28 additions & 18 deletions source/src_pw/klist.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -80,12 +80,16 @@ void K_Vectors::set(
ModuleBase::WARNING_QUIT("K_Vectors::set","Something wrong while reading KPOINTS.");
}

// output kpoints file
std::string skpt1="";
std::string skpt2="";

// (2)
//only berry phase need all kpoints including time-reversal symmetry!
//if symm_flag is not set, only time-reversal symmetry would be considered.
if(!berryphase::berry_phase_flag && ModuleSymmetry::Symmetry::symm_flag != -1)
{
this->ibz_kpoint(symm, ModuleSymmetry::Symmetry::symm_flag);
this->ibz_kpoint(symm, ModuleSymmetry::Symmetry::symm_flag,skpt1);
if(ModuleSymmetry::Symmetry::symm_flag || is_mp)
{
this->update_use_ibz();
Expand All @@ -94,7 +98,17 @@ void K_Vectors::set(
}

// (3)
this->set_both_kvec(reciprocal_vec, latvec);
this->set_both_kvec(reciprocal_vec, latvec, skpt2);

if(GlobalV::MY_RANK==0)
{
// output kpoints file
std::stringstream skpt;
skpt << GlobalV::global_readin_dir << "kpoints" ;
std::ofstream ofkpt( skpt.str().c_str()); // clear kpoints
ofkpt << skpt2 <<skpt1;
ofkpt.close();
}

int deg = 0;
if(GlobalV::NSPIN == 1)
Expand Down Expand Up @@ -534,7 +548,7 @@ void K_Vectors::update_use_ibz( void )
return;
}

void K_Vectors::ibz_kpoint(const ModuleSymmetry::Symmetry &symm, bool use_symm)
void K_Vectors::ibz_kpoint(const ModuleSymmetry::Symmetry &symm, bool use_symm,std::string& skpt)
{
if (GlobalV::MY_RANK!=0) return;
ModuleBase::TITLE("K_Vectors", "ibz_kpoint");
Expand Down Expand Up @@ -602,11 +616,9 @@ void K_Vectors::ibz_kpoint(const ModuleSymmetry::Symmetry &symm, bool use_symm)
// }
// output in kpoints file
std::stringstream ss;
ss << GlobalV::global_readin_dir << "kpoints" ;
std::ofstream ofkpt( ss.str().c_str()); // clear kpoints
ofkpt<< " " << std::setw(40) <<"nkstot" << " = " << nkstot
ss << " " << std::setw(40) <<"nkstot" << " = " << nkstot
<< std::setw(66) << "ibzkpt" << std::endl;
ofkpt << " " << std::setw(8) << "KPT"
ss << " " << std::setw(8) << "KPT"
<< std::setw(20) << "DirectX"
<< std::setw(20) << "DirectY"
<< std::setw(20) << "DirectZ"
Expand All @@ -619,7 +631,7 @@ void K_Vectors::ibz_kpoint(const ModuleSymmetry::Symmetry &symm, bool use_symm)
for (int i = 0; i < nkstot; ++i)
{
// output in kpoints file
ofkpt << " "
ss << " "
<< std::setw(8) << i+1
<< std::setw(20) << this->kvec_d[i].x
<< std::setw(20) << this->kvec_d[i].y
Expand Down Expand Up @@ -660,7 +672,7 @@ void K_Vectors::ibz_kpoint(const ModuleSymmetry::Symmetry &symm, bool use_symm)
// but is already in the ibz_kpoint list.
// so the weight need to +1;
// output in kpoints file
ofkpt << std::setw(8) << k+1
ss << std::setw(8) << k+1
<< std::setw(20) << this->kvec_d_ibz[k].x
<< std::setw(20) << this->kvec_d_ibz[k].y
<< std::setw(20) << this->kvec_d_ibz[k].z << std::endl;
Expand All @@ -679,7 +691,7 @@ void K_Vectors::ibz_kpoint(const ModuleSymmetry::Symmetry &symm, bool use_symm)
//nkstot_ibz indicate the index of ibz kpoint.
this->kvec_d_ibz[nkstot_ibz] = kvec_rot;
// output in kpoints file
ofkpt << std::setw(8) << nkstot_ibz+1
ss << std::setw(8) << nkstot_ibz+1
<< std::setw(20) << this->kvec_d_ibz[nkstot_ibz].x
<< std::setw(20) << this->kvec_d_ibz[nkstot_ibz].y
<< std::setw(20) << this->kvec_d_ibz[nkstot_ibz].z << std::endl;
Expand Down Expand Up @@ -722,7 +734,7 @@ void K_Vectors::ibz_kpoint(const ModuleSymmetry::Symmetry &symm, bool use_symm)
// BLOCK_HERE("check k point");
}

ofkpt.close();
skpt = ss.str();

ModuleBase::GlobalFunc::OUT(GlobalV::ofs_running,"nkstot_ibz",nkstot_ibz);

Expand All @@ -744,7 +756,7 @@ void K_Vectors::ibz_kpoint(const ModuleSymmetry::Symmetry &symm, bool use_symm)
}


void K_Vectors::set_both_kvec(const ModuleBase::Matrix3 &G, const ModuleBase::Matrix3 &R)
void K_Vectors::set_both_kvec(const ModuleBase::Matrix3 &G, const ModuleBase::Matrix3 &R,std::string& skpt)
{

if(GlobalV::FINAL_SCF) //LiuXh add 20180606
Expand Down Expand Up @@ -830,22 +842,20 @@ void K_Vectors::set_both_kvec(const ModuleBase::Matrix3 &G, const ModuleBase::Ma
if(GlobalV::MY_RANK==0)
{
std::stringstream ss;
ss << GlobalV::global_readin_dir << "kpoints" ;
std::ofstream ofkpt( ss.str().c_str(), ios::app ); // clear kpoints
ofkpt << " " << std::setw(40) <<"nkstot now" << " = " << nkstot << std::endl;
ofkpt << " " << std::setw(8) << "KPT" << std::setw(20) << "DirectX"
ss << " " << std::setw(40) <<"nkstot now" << " = " << nkstot << std::endl;
ss << " " << std::setw(8) << "KPT" << std::setw(20) << "DirectX"
<< std::setw(20) << "DirectY" << std::setw(20) << "DirectZ"
<< std::setw(20) << "Weight" << std::endl;
for (int ik=0; ik<nkstot; ik++)
{
ofkpt << " " << std::setw(8) << ik+1
ss << " " << std::setw(8) << ik+1
<< std::setw(20) << this->kvec_d[ik].x
<< std::setw(20) << this->kvec_d[ik].y
<< std::setw(20) << this->kvec_d[ik].z
<< std::setw(20) << this->wk[ik] << std::endl;
// << std::setw(10) << this->ibz2bz[ik] << std::endl;
}
ofkpt.close();
skpt = ss.str();
}

return;
Expand Down
4 changes: 2 additions & 2 deletions source/src_pw/klist.h
Original file line number Diff line number Diff line change
Expand Up @@ -38,7 +38,7 @@ class K_Vectors
const ModuleBase::Matrix3 &reciprocal_vec,
const ModuleBase::Matrix3 &latvec);

void ibz_kpoint( const ModuleSymmetry::Symmetry &symm, bool use_symm);
void ibz_kpoint( const ModuleSymmetry::Symmetry &symm, bool use_symm,std::string& skpt);
//LiuXh add 20180515
void set_after_vc(
const ModuleSymmetry::Symmetry &symm,
Expand Down Expand Up @@ -66,7 +66,7 @@ class K_Vectors

// step 2 : set both kvec and kved; normalize weight
void update_use_ibz( void );
void set_both_kvec(const ModuleBase::Matrix3 &G,const ModuleBase::Matrix3 &R);
void set_both_kvec(const ModuleBase::Matrix3 &G,const ModuleBase::Matrix3 &R, std::string& skpt);
void normalize_wk( const int &degspin );

// step 3 : mpi kpoints information.
Expand Down

0 comments on commit f145042

Please sign in to comment.