|
ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
|
#include <sto_iter.h>
Public Member Functions | |
| Stochastic_Iter () | |
| ~Stochastic_Iter () | |
| void | init (K_Vectors *pkv_in, ModulePW::PW_Basis_K *wfc_basis, Stochastic_WF< T, Device > &stowf, StoChe< Real, Device > &stoche, hamilt::HamiltSdftPW< T, Device > *p_hamilt_sto) |
| init for iteration process of SDFT | |
| void | sum_stoeband (Stochastic_WF< T, Device > &stowf, elecstate::ElecStatePW< T, Device > *pes, hamilt::Hamilt< T, Device > *pHamilt, ModulePW::PW_Basis_K *wfc_basis) |
| sum demet and eband energies for each k point and each band | |
| void | cal_storho (const UnitCell &ucell, Stochastic_WF< T, Device > &stowf, elecstate::ElecStatePW< T, Device > *pes, ModulePW::PW_Basis_K *wfc_basis) |
| calculate the density | |
| double | calne (elecstate::ElecState *pes) |
| calculate total number of electrons | |
| void | itermu (const int iter, elecstate::ElecState *pes) |
| solve ne(mu) = ne_target and get chemical potential mu | |
| void | orthog (const int &ik, psi::Psi< T, Device > &psi, Stochastic_WF< T, Device > &stowf) |
| orthogonalize stochastic wave functions with KS wave functions | |
| void | checkemm (const int &ik, const int istep, const int iter, Stochastic_WF< T, Device > &stowf) |
| check emax and emin | |
| void | check_precision (const double ref, const double thr, const std::string info) |
| check precision of Chebyshev expansion | |
| void | calHsqrtchi (Stochastic_WF< T, Device > &stowf) |
| void | calPn (const int &ik, Stochastic_WF< T, Device > &stowf) |
| void | calTnchi_ik (const int &ik, Stochastic_WF< T, Device > &stowf) |
Public Attributes | |
| ModuleBase::Chebyshev< double, Device > * | p_che = nullptr |
| Sto_Func< double > | stofunc |
| hamilt::HamiltSdftPW< T, Device > * | p_hamilt_sto = nullptr |
| double | mu0 |
| bool | change |
| double | targetne =0.0 |
| Real * | spolyv = nullptr |
| Real * | spolyv_cpu = nullptr |
| int * | nchip = nullptr |
| bool | check = false |
| double | th_ne =0.0 |
| double | KS_ne =0.0 |
| int | method |
Private Member Functions | |
| void | dot (const int &n, const Real *x, const int &incx, const Real *y, const int &incy, Real &result) |
| return cpu dot result | |
Private Attributes | |
| K_Vectors * | pkv =nullptr |
| const Device * | ctx = {} |
| const base_device::DEVICE_CPU * | cpu_ctx = {} |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
| Stochastic_Iter< T, Device >::Stochastic_Iter | ( | ) |
| Stochastic_Iter< T, Device >::~Stochastic_Iter | ( | ) |
| void Stochastic_Iter< T, Device >::cal_storho | ( | const UnitCell & | ucell, |
| Stochastic_WF< T, Device > & | stowf, | ||
| elecstate::ElecStatePW< T, Device > * | pes, | ||
| ModulePW::PW_Basis_K * | wfc_basis | ||
| ) |
calculate the density
| ucell | reference to unit cell |
| stowf | stochastic wave function |
| pes | elecstate |
| wfc_basis | wfc pw basis |
| void Stochastic_Iter< T, Device >::calHsqrtchi | ( | Stochastic_WF< T, Device > & | stowf | ) |
| double Stochastic_Iter< T, Device >::calne | ( | elecstate::ElecState * | pes | ) |
calculate total number of electrons
| pes | elecstate |
| void Stochastic_Iter< T, Device >::calPn | ( | const int & | ik, |
| Stochastic_WF< T, Device > & | stowf | ||
| ) |
| void Stochastic_Iter< T, Device >::calTnchi_ik | ( | const int & | ik, |
| Stochastic_WF< T, Device > & | stowf | ||
| ) |
| void Stochastic_Iter< T, Device >::check_precision | ( | const double | ref, |
| const double | thr, | ||
| const std::string | info | ||
| ) |
check precision of Chebyshev expansion
| ref | reference value |
| thr | threshold |
| info | information |
| void Stochastic_Iter< T, Device >::checkemm | ( | const int & | ik, |
| const int | istep, | ||
| const int | iter, | ||
| Stochastic_WF< T, Device > & | stowf | ||
| ) |
check emax and emin
| ik | k point index |
| istep | ion step index |
| iter | scf iteration index |
| stowf | stochastic wave functions |
|
private |
return cpu dot result
| x | [Device] |
| y | [Device] |
| result | [CPU] dot result |
| void Stochastic_Iter< T, Device >::init | ( | K_Vectors * | pkv_in, |
| ModulePW::PW_Basis_K * | wfc_basis, | ||
| Stochastic_WF< T, Device > & | stowf, | ||
| StoChe< Real, Device > & | stoche, | ||
| hamilt::HamiltSdftPW< T, Device > * | p_hamilt_sto | ||
| ) |
init for iteration process of SDFT
| method_in | 1: slow 2: fast but cost much memories |
| pkv_in | K_Vectors |
| wfc_basis | wfc pw basis |
| stowf | stochastic wave function |
| stoche | Chebyshev expansion for sDFT |
| p_hamilt_sto | hamiltonian for sDFT |
| void Stochastic_Iter< T, Device >::itermu | ( | const int | iter, |
| elecstate::ElecState * | pes | ||
| ) |
solve ne(mu) = ne_target and get chemical potential mu
| iter | scf iteration index |
| pes | elecstate |
| void Stochastic_Iter< T, Device >::orthog | ( | const int & | ik, |
| psi::Psi< T, Device > & | psi, | ||
| Stochastic_WF< T, Device > & | stowf | ||
| ) |
orthogonalize stochastic wave functions with KS wave functions
| ik | k point index |
| psi | KS wave functions |
| stowf | stochastic wave functions |
| void Stochastic_Iter< T, Device >::sum_stoeband | ( | Stochastic_WF< T, Device > & | stowf, |
| elecstate::ElecStatePW< T, Device > * | pes, | ||
| hamilt::Hamilt< T, Device > * | pHamilt, | ||
| ModulePW::PW_Basis_K * | wfc_basis | ||
| ) |
sum demet and eband energies for each k point and each band
| stowf | stochastic wave function |
| pes | elecstate |
| pHamilt | hamiltonian |
| wfc_basis | wfc pw basis |
| bool Stochastic_Iter< T, Device >::change |
| bool Stochastic_Iter< T, Device >::check = false |
|
private |
|
private |
| double Stochastic_Iter< T, Device >::KS_ne =0.0 |
| int Stochastic_Iter< T, Device >::method |
| double Stochastic_Iter< T, Device >::mu0 |
| int* Stochastic_Iter< T, Device >::nchip = nullptr |
| ModuleBase::Chebyshev<double, Device>* Stochastic_Iter< T, Device >::p_che = nullptr |
| hamilt::HamiltSdftPW<T, Device>* Stochastic_Iter< T, Device >::p_hamilt_sto = nullptr |
|
private |
| Real* Stochastic_Iter< T, Device >::spolyv = nullptr |
| Real* Stochastic_Iter< T, Device >::spolyv_cpu = nullptr |
| Sto_Func<double> Stochastic_Iter< T, Device >::stofunc |
| double Stochastic_Iter< T, Device >::targetne =0.0 |
| double Stochastic_Iter< T, Device >::th_ne =0.0 |