Skip to content

Commit 07dc729

Browse files
committed
Fix: using psi rather than evc in winput::out_spillage==1 or 2
1 parent e8d5d5b commit 07dc729

File tree

3 files changed

+11
-9
lines changed

3 files changed

+11
-9
lines changed

source/run_pw.cpp

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -96,7 +96,7 @@ void Run_pw::plane_wave_line(ModuleESolver::ESolver *p_esolver)
9696
if ( winput::out_spillage <= 2 )
9797
{
9898
Numerical_Basis numerical_basis;
99-
numerical_basis.output_overlap(GlobalC::wf.evc);
99+
numerical_basis.output_overlap(GlobalC::wf.psi[0]);
100100
ModuleBase::GlobalFunc::DONE(GlobalV::ofs_running,"BASIS OVERLAP (Q and S) GENERATION.");
101101
}
102102
}

source/src_io/numerical_basis.cpp

Lines changed: 6 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -43,7 +43,7 @@ void Numerical_Basis::start_from_file_k( const int &ik, ModuleBase::ComplexMatri
4343
}
4444

4545
// The function is called in run_fp.cpp.
46-
void Numerical_Basis::output_overlap( const ModuleBase::ComplexMatrix *psi)
46+
void Numerical_Basis::output_overlap( const psi::Psi<std::complex<double>> &psi)
4747
{
4848
ModuleBase::TITLE("Numerical_Basis","output_overlap");
4949
ModuleBase::GlobalFunc::NEW_PART("Overlap Data For Spillage Minimization");
@@ -97,7 +97,8 @@ void Numerical_Basis::output_overlap( const ModuleBase::ComplexMatrix *psi)
9797
GlobalV::ofs_running << " --------------------------------------------------------" << std::endl;
9898

9999
// search for all k-points.
100-
overlap_Q[ik] = this->cal_overlap_Q(ik, npw, psi[ik], static_cast<double>(derivative_order));
100+
psi.fix_k(ik);
101+
overlap_Q[ik] = this->cal_overlap_Q(ik, npw, psi, static_cast<double>(derivative_order));
101102
ModuleBase::GlobalFunc::DONE(GlobalV::ofs_running,"cal_overlap_Q");
102103

103104
// (2) generate Sq matrix if necessary.
@@ -140,7 +141,7 @@ void Numerical_Basis::output_overlap( const ModuleBase::ComplexMatrix *psi)
140141
ModuleBase::ComplexArray Numerical_Basis::cal_overlap_Q(
141142
const int &ik,
142143
const int &np,
143-
const ModuleBase::ComplexMatrix &psi,
144+
const psi::Psi<std::complex<double>> &psi,
144145
const double derivative_order) const
145146
{
146147
ModuleBase::TITLE("Numerical_Basis","cal_overlap_Q");
@@ -335,7 +336,7 @@ ModuleBase::ComplexArray Numerical_Basis::cal_overlap_Sq(
335336

336337
// Peize Lin add for dpsi 2020.04.23
337338
ModuleBase::matrix Numerical_Basis::cal_overlap_V(
338-
const ModuleBase::ComplexMatrix *psi,
339+
const psi::Psi<std::complex<double>> &psi,
339340
const double derivative_order)
340341
{
341342
ModuleBase::matrix overlap_V(GlobalC::kv.nks, GlobalV::NBANDS);
@@ -349,7 +350,7 @@ ModuleBase::matrix Numerical_Basis::cal_overlap_V(
349350

350351
for(int ib=0; ib<GlobalV::NBANDS; ++ib)
351352
for(int ig=0; ig<GlobalC::kv.ngk[ik]; ++ig)
352-
overlap_V(ik,ib)+= norm(psi[ik](ib,ig)) * gpow[ig];
353+
overlap_V(ik,ib)+= norm(psi(ik,ib,ig)) * gpow[ig];
353354
}
354355
return overlap_V;
355356
}

source/src_io/numerical_basis.h

Lines changed: 4 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -14,6 +14,7 @@
1414
#include "../module_base/complexmatrix.h"
1515
#include "bessel_basis.h"
1616
#include <vector>
17+
#include "module_psi/psi.h"
1718
//==========================================================
1819
// CLASS :
1920
// NAME : Numerical_Basis
@@ -25,7 +26,7 @@ class Numerical_Basis
2526
~Numerical_Basis();
2627

2728
void start_from_file_k( const int &ik, ModuleBase::ComplexMatrix &psi);
28-
void output_overlap( const ModuleBase::ComplexMatrix *psi);
29+
void output_overlap( const psi::Psi<std::complex<double>> &psi);
2930

3031
private:
3132

@@ -39,15 +40,15 @@ class Numerical_Basis
3940
void numerical_atomic_wfc(const int &ik,const int &np,ModuleBase::ComplexMatrix &psi);
4041

4142
ModuleBase::ComplexArray cal_overlap_Q(
42-
const int &ik, const int &np, const ModuleBase::ComplexMatrix &psi,
43+
const int &ik, const int &np, const psi::Psi<std::complex<double>> &psi,
4344
const double derivative_order ) const;
4445

4546
ModuleBase::ComplexArray cal_overlap_Sq(
4647
const int &ik,
4748
const int &np,
4849
const double derivative_order ) const;
4950

50-
static ModuleBase::matrix cal_overlap_V(const ModuleBase::ComplexMatrix *psi, const double derivative_order);
51+
static ModuleBase::matrix cal_overlap_V(const psi::Psi<std::complex<double>> &psi, const double derivative_order);
5152

5253
ModuleBase::realArray cal_flq(const int ik, const std::vector<ModuleBase::Vector3<double>> &gk) const;
5354

0 commit comments

Comments
 (0)