SIRIUS 7.5.0
Electronic structure library and applications
Public Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
sirius::K_point< T > Class Template Reference

K-point related variables and methods. More...

#include <k_point.hpp>

Public Member Functions

 K_point (Simulation_context &ctx__, r3::vector< double > vk__, double weight__)
 Constructor. More...
 
 K_point (Simulation_context &ctx__, std::shared_ptr< fft::Gvec > gkvec__, double weight__)
 Constructor. More...
 
void initialize ()
 Initialize the k-point related arrays and data. More...
 
void update ()
 Update the reciprocal lattice vectors of the G+k array. More...
 
void generate_fv_states ()
 Generate first-variational states from eigen-vectors. More...
 
void generate_spinor_wave_functions ()
 Generate two-component spinor wave functions. More...
 
void generate_atomic_wave_functions (std::vector< int > atoms__, std::function< basis_functions_index const *(int)> indexb__, Radial_integrals_atomic_wf< false > const &ri__, wf::Wave_functions< T > &wf__)
 Generate plane-wave coefficients of the atomic wave-functions. More...
 
void generate_hubbard_orbitals ()
 
void save (std::string const &name__, int id__) const
 Save data to HDF5 file. More...
 
void load (HDF5_tree h5in, int id)
 
void get_fv_eigen_vectors (sddk::mdarray< std::complex< T >, 2 > &fv_evec__) const
 Collect distributed first-variational vectors into a global array. More...
 
void get_sv_eigen_vectors (sddk::mdarray< std::complex< T >, 2 > &sv_evec__) const
 Collect distributed second-variational vectors into a global array. More...
 
auto const & gkvec () const
 
auto gkvec_sptr () const
 Return shared pointer to gkvec object. More...
 
int num_gkvec () const
 Total number of G+k vectors within the cutoff distance. More...
 
int num_gkvec_loc () const
 Local number of G+k vectors in case of flat distribution. More...
 
int num_occupied_bands (int ispn__=-1) const
 Get the number of occupied bands for each spin channel. More...
 
double band_energy (int j__, int ispn__) const
 Get band energy. More...
 
void band_energy (int j__, int ispn__, double e__)
 Set band energy. More...
 
auto band_energies (int ispn__) const
 
double band_occupancy (int j__, int ispn__) const
 Get band occupancy. More...
 
void band_occupancy (int j__, int ispn__, double occ__)
 Set band occupancy. More...
 
double fv_eigen_value (int i) const
 
void set_fv_eigen_values (double *eval__)
 
double weight () const
 Return weight of k-point. More...
 
auto & fv_states ()
 
wf::Wave_functions< T > const & spinor_wave_functions () const
 
wf::Wave_functions< T > & spinor_wave_functions ()
 
auto & spinor_wave_functions2 ()
 
auto const & atomic_wave_functions_S () const
 
auto & atomic_wave_functions_S ()
 
auto const & atomic_wave_functions () const
 Return the initial atomic orbitals used to compute the hubbard wave functions. More...
 
auto & atomic_wave_functions ()
 Return the initial atomic orbitals used to compute the hubbard wave functions. More...
 
auto const & hubbard_wave_functions_S () const
 Return the actual hubbard wave functions used in the calculations. More...
 
auto & singular_components ()
 
auto vk () const
 
int gklo_basis_size () const
 Basis size of LAPW+lo method. More...
 
int num_gkvec_row () const
 Local number of G+k vectors for each MPI rank in the row of the 2D MPI grid. More...
 
int num_lo_row () const
 Local number of local orbitals for each MPI rank in the row of the 2D MPI grid. More...
 
int gklo_basis_size_row () const
 Local number of basis functions for each MPI rank in the row of the 2D MPI grid. More...
 
int num_gkvec_col () const
 Local number of G+k vectors for each MPI rank in the column of the 2D MPI grid. More...
 
int num_lo_col () const
 Local number of local orbitals for each MPI rank in the column of the 2D MPI grid. More...
 
int gklo_basis_size_col () const
 Local number of basis functions for each MPI rank in the column of the 2D MPI grid. More...
 
lo_basis_descriptor const & lo_basis_descriptor_col (int idx) const
 
lo_basis_descriptor const & lo_basis_descriptor_row (int idx) const
 
int num_ranks_row () const
 
int rank_row () const
 
int num_ranks_col () const
 
int rank_col () const
 
int num_ranks () const
 Number of MPI ranks for a given k-point. More...
 
int rank () const
 
int num_atom_lo_cols (int ia) const
 Return number of lo columns for a given atom. More...
 
int lo_col (int ia, int i) const
 Return local index (for the current MPI rank) of a column for a given atom and column index within an atom. More...
 
int num_atom_lo_rows (int ia) const
 Return number of lo rows for a given atom. More...
 
int lo_row (int ia, int i) const
 Return local index (for the current MPI rank) of a row for a given atom and row index within an atom. More...
 
auto & fv_eigen_vectors ()
 
auto & fv_eigen_vectors_slab ()
 
auto & sv_eigen_vectors (int ispn)
 
auto & fd_eigen_vectors ()
 
void bypass_sv ()
 
auto const & alm_coeffs_row () const
 
auto const & alm_coeffs_col () const
 
auto const & alm_coeffs_loc () const
 
auto const & comm () const
 
auto const & comm_row () const
 
auto const & comm_col () const
 
auto beta_projectors () -> Beta_projectors< T > &
 
auto beta_projectors () const -> const Beta_projectors< T > &
 
auto & beta_projectors_row ()
 
auto & beta_projectors_col ()
 
auto const & ctx () const
 
std::ostream & out (int level__) const
 Return stdout stream for high verbosity output. More...
 
auto & spfft_transform () const
 
auto const & gkvec_fft () const
 
auto gkvec_fft_sptr () const
 
auto const & gkvec_col () const
 
auto const & gkvec_row () const
 

Private Member Functions

void generate_gklo_basis ()
 Generate G+k and local orbital basis sets. More...
 
void generate_gkvec (double gk_cutoff__)
 Find G+k vectors within the cutoff. More...
 
int get_ispn (int ispn__) const
 
void init0 ()
 

Private Attributes

Simulation_contextctx_
 Simulation context. More...
 
Unit_cell const & unit_cell_
 Unit cell object. More...
 
r3::vector< double > vk_
 Fractional k-point coordinates. More...
 
double weight_ {1.0}
 Weight of k-point. More...
 
mpi::Communicator const & comm_
 Communicator for parallelization inside k-point. More...
 
std::shared_ptr< fft::Gvecgkvec_
 List of G-vectors with |G+k| < cutoff. More...
 
std::shared_ptr< fft::Gvecgkvec_row_
 
std::shared_ptr< fft::Gvecgkvec_col_
 
std::shared_ptr< fft::Gvec_fftgkvec_partition_
 G-vector distribution for the FFT transformation. More...
 
std::unique_ptr< fft::spfft_transform_type< T > > spfft_transform_
 
sddk::mdarray< double, 1 > fv_eigen_values_
 First-variational eigen values. More...
 
la::dmatrix< std::complex< T > > fv_eigen_vectors_
 First-variational eigen vectors, distributed over 2D BLACS grid. More...
 
std::unique_ptr< wf::Wave_functions< T > > fv_eigen_vectors_slab_
 First-variational eigen vectors, distributed in slabs. More...
 
std::unique_ptr< wf::Wave_functions< T > > singular_components_
 Lowest eigen-vectors of the LAPW overlap matrix with small aigen-values. More...
 
std::array< la::dmatrix< std::complex< T > >, 2 > sv_eigen_vectors_
 Second-variational eigen vectors. More...
 
sddk::mdarray< std::complex< T >, 2 > fd_eigen_vectors_
 Full-diagonalization eigen vectors. More...
 
std::unique_ptr< wf::Wave_functions< T > > fv_states_ {nullptr}
 First-variational states. More...
 
std::unique_ptr< wf::Wave_functions< T > > spinor_wave_functions_ {nullptr}
 Two-component (spinor) wave functions describing the bands. More...
 
std::unique_ptr< wf::Wave_functions< T > > atomic_wave_functions_ {nullptr}
 Pseudopotential atmoic wave-functions (not orthogonalized). More...
 
std::unique_ptr< wf::Wave_functions< T > > atomic_wave_functions_S_ {nullptr}
 Pseudopotential atmoic wave-functions (not orthogonalized) with S-operator applied. More...
 
std::unique_ptr< wf::Wave_functions< T > > hubbard_wave_functions_ {nullptr}
 Hubbard wave functions. More...
 
std::unique_ptr< wf::Wave_functions< T > > hubbard_wave_functions_S_ {nullptr}
 Hubbard wave functions with S-operator applied. More...
 
sddk::mdarray< double, 2 > band_occupancies_
 Band occupation numbers. More...
 
sddk::mdarray< double, 2 > band_energies_
 Band energies. More...
 
std::unique_ptr< Matching_coefficientsalm_coeffs_row_ {nullptr}
 LAPW matching coefficients for the row G+k vectors. More...
 
std::unique_ptr< Matching_coefficientsalm_coeffs_col_ {nullptr}
 LAPW matching coefficients for the column G+k vectors. More...
 
std::unique_ptr< Matching_coefficientsalm_coeffs_loc_ {nullptr}
 LAPW matching coefficients for the local set G+k vectors. More...
 
int num_gkvec_row_ {0}
 Number of G+k vectors distributed along rows of MPI grid. More...
 
int num_gkvec_col_ {0}
 Number of G+k vectors distributed along columns of MPI grid. More...
 
std::vector< lo_basis_descriptorlo_basis_descriptors_row_
 Basis descriptors distributed between rows of the 2D MPI grid. More...
 
std::vector< lo_basis_descriptorlo_basis_descriptors_col_
 Basis descriptors distributed between columns of the 2D MPI grid. More...
 
std::vector< std::vector< int > > atom_lo_cols_
 List of columns of the Hamiltonian and overlap matrix lo block (local index) for a given atom. More...
 
std::vector< std::vector< int > > atom_lo_rows_
 list of rows of the Hamiltonian and overlap matrix lo block (local index) for a given atom More...
 
std::vector< std::complex< double > > zil_
 Imaginary unit to the power of l. More...
 
std::vector< int > l_by_lm_
 Mapping between lm and l. More...
 
int rank_col_
 Column rank of the processors of ScaLAPACK/ELPA diagonalization grid. More...
 
int num_ranks_col_
 Number of processors along the columns of the diagonalization grid. More...
 
int rank_row_
 Row rank of the processors of ScaLAPACK/ELPA diagonalization grid. More...
 
int num_ranks_row_
 Number of processors along the rows of the diagonalization grid. More...
 
std::unique_ptr< Beta_projectors< T > > beta_projectors_ {nullptr}
 Beta projectors for a local set of G+k vectors. More...
 
std::unique_ptr< Beta_projectors< T > > beta_projectors_row_ {nullptr}
 Beta projectors for row G+k vectors. More...
 
std::unique_ptr< Beta_projectors< T > > beta_projectors_col_ {nullptr}
 Beta projectors for column G+k vectors. More...
 
mpi::Communicator const & comm_row_
 Communicator between(!!) rows. More...
 
mpi::Communicator const & comm_col_
 Communicator between(!!) columns. More...
 
std::array< int, 2 > ispn_map_ {0, -1}
 

Friends

class K_point_set
 

Detailed Description

template<typename T>
class sirius::K_point< T >

K-point related variables and methods.

Wave-function storage
First-variational eigen vectors
Template Parameters
TPrecision of the wave-functions (float or double).

Definition at line 43 of file k_point.hpp.

Constructor & Destructor Documentation

◆ K_point() [1/2]

template<typename T >
sirius::K_point< T >::K_point ( Simulation_context ctx__,
r3::vector< double >  vk__,
double  weight__ 
)
inline

Constructor.

Definition at line 219 of file k_point.hpp.

◆ K_point() [2/2]

template<typename T >
sirius::K_point< T >::K_point ( Simulation_context ctx__,
std::shared_ptr< fft::Gvec gkvec__,
double  weight__ 
)
inline

Constructor.

Definition at line 238 of file k_point.hpp.

Member Function Documentation

◆ generate_gklo_basis()

template<typename T >
void sirius::K_point< T >::generate_gklo_basis
private

Generate G+k and local orbital basis sets.

Definition at line 763 of file k_point.cpp.

◆ generate_gkvec()

template<typename T >
void sirius::K_point< T >::generate_gkvec ( double  gk_cutoff__)
private

Find G+k vectors within the cutoff.

Definition at line 340 of file k_point.cpp.

◆ get_ispn()

template<typename T >
int sirius::K_point< T >::get_ispn ( int  ispn__) const
inlineprivate

Definition at line 193 of file k_point.hpp.

◆ init0()

template<typename T >
void sirius::K_point< T >::init0 ( )
inlineprivate

Definition at line 201 of file k_point.hpp.

◆ initialize()

template<typename T >
void sirius::K_point< T >::initialize

Initialize the k-point related arrays and data.

Definition at line 33 of file k_point.cpp.

◆ update()

template<typename T >
void sirius::K_point< T >::update

Update the reciprocal lattice vectors of the G+k array.

Definition at line 393 of file k_point.cpp.

◆ generate_fv_states()

template<typename T >
template void sirius::K_point< T >::generate_fv_states ( )

Generate first-variational states from eigen-vectors.

APW+lo basis \( \varphi_{\mu {\bf k}}({\bf r}) = \{ \varphi_{\bf G+k}({\bf r}), \varphi_{j{\bf k}}({\bf r}) \} \) is used to expand first-variational wave-functions:

\[ \psi_{i{\bf k}}({\bf r}) = \sum_{\mu} c_{\mu i}^{\bf k} \varphi_{\mu \bf k}({\bf r}) = \sum_{{\bf G}}c_{{\bf G} i}^{\bf k} \varphi_{\bf G+k}({\bf r}) + \sum_{j}c_{j i}^{\bf k}\varphi_{j{\bf k}}({\bf r}) \]

Inside muffin-tins the expansion is converted into the following form:

\[ \psi_{i {\bf k}}({\bf r})= \begin{array}{ll} \displaystyle \sum_{L} \sum_{\lambda=1}^{N_{\ell}^{\alpha}} F_{L \lambda}^{i {\bf k},\alpha}f_{\ell \lambda}^{\alpha}(r) Y_{\ell m}(\hat {\bf r}) & {\bf r} \in MT_{\alpha} \end{array} \]

Thus, the total number of coefficients representing a wave-funstion is equal to the number of muffin-tin basis functions of the form \( f_{\ell \lambda}^{\alpha}(r) Y_{\ell m}(\hat {\bf r}) \) plust the number of G+k plane waves. First-variational states are obtained from the first-variational eigen-vectors and LAPW matching coefficients.

APW part:

\[ \psi_{\xi j}^{\bf k} = \sum_{{\bf G}} Z_{{\bf G} j}^{\bf k} * A_{\xi}({\bf G+k}) \]

Definition at line 31 of file generate_fv_states.cpp.

◆ generate_spinor_wave_functions()

template<typename T >
template void sirius::K_point< T >::generate_spinor_wave_functions ( )

Generate two-component spinor wave functions.

In case of second-variational diagonalization spinor wave-functions are generated from the first-variational states and second-variational eigen-vectors.

Definition at line 30 of file generate_spinor_wave_functions.cpp.

◆ generate_atomic_wave_functions()

template<typename T >
void sirius::K_point< T >::generate_atomic_wave_functions ( std::vector< int >  atoms__,
std::function< basis_functions_index const *(int)>  indexb__,
Radial_integrals_atomic_wf< false > const &  ri__,
wf::Wave_functions< T > &  wf__ 
)

Generate plane-wave coefficients of the atomic wave-functions.

Plane-wave coefficients of the atom-centered wave-functions \( \varphi^{\alpha}_{\ell m}({\bf r}) = \varphi^{\alpha}_{\ell}(r)R_{\ell m}(\theta, \phi) \) are computed in the following way:

\[ \varphi^{\alpha}_{\ell m}({\bf q}) = \frac{1}{\sqrt{\Omega}} \int e^{-i{\bf q}{\bf r}} \varphi^{\alpha}_{\ell m}({\bf r} - {\bf r}_{\alpha}) d{\bf r} = \frac{e^{-i{\bf q}{\bf r}_{\alpha}}}{\sqrt{\Omega}} \int e^{-i{\bf q}{\bf r}} \varphi^{\alpha}_{\ell}(r)R_{\ell m}(\theta, \phi) r^2 \sin\theta dr d\theta d\phi \]

where \( {\bf q} = {\bf G+k} \). Using the expansion of the plane wave in terms of spherical Bessel functions and real spherical harmonics:

\[ e^{-i{\bf q}{\bf r}}=4\pi \sum_{\ell m} (-i)^\ell j_{\ell}(q r)R_{\ell m}({\bf \hat q})R_{\ell m}({\bf \hat r}) \]

we arrive to the following expression:

\[ \varphi^{\alpha}_{\ell m}({\bf q}) = e^{-i{\bf q}{\bf r}_{\alpha}} \frac{4\pi}{\sqrt{\Omega}} (-i)^\ell R_{\ell m}({\bf q}) \int \varphi^{\alpha}_{\ell}(r) j_{\ell}(q r) r^2 dr \]

Note
In the current implementation wave-functions are generated as scalars (without spin index). Spinor atomic wave-functions might be necessary in future for the more advanced LDA+U implementation.
Parameters
[in]atomsList of atoms, for which the wave-functions are generated.
[in]indexbLambda function that returns index of the basis functions for each atom type.
[in]riRadial integrals of the product of sperical Bessel functions and atomic functions.
[out]wfResulting wave-functions for the list of atoms. Output wave-functions must have sufficient storage space.

Definition at line 678 of file k_point.cpp.

◆ generate_hubbard_orbitals()

template<typename T >
void sirius::K_point< T >::generate_hubbard_orbitals

Definition at line 189 of file k_point.cpp.

◆ save()

template<typename T >
void sirius::K_point< T >::save ( std::string const &  name__,
int  id__ 
) const

Save data to HDF5 file.

The following HDF5 data structure is created:

/K_point_set/ik/vk
/K_point_set/ik/band_energies
/K_point_set/ik/band_occupancies
/K_point_set/ik/gkvec
/K_point_set/ik/gvec
/K_point_set/ik/bands/ibnd/spinor_wave_function/ispn/pw
/K_point_set/ik/bands/ibnd/spinor_wave_function/ispn/mt
  • store wave-functions *‍/

Definition at line 539 of file k_point.cpp.

◆ load()

template<typename T >
void sirius::K_point< T >::load ( HDF5_tree  h5in,
int  id 
)

Definition at line 605 of file k_point.cpp.

◆ get_fv_eigen_vectors()

template<typename T >
void sirius::K_point< T >::get_fv_eigen_vectors ( sddk::mdarray< std::complex< T >, 2 > &  fv_evec__) const

Collect distributed first-variational vectors into a global array.

Definition at line 425 of file k_point.cpp.

◆ get_sv_eigen_vectors()

template<typename T >
void sirius::K_point< T >::get_sv_eigen_vectors ( sddk::mdarray< std::complex< T >, 2 > &  sv_evec__) const
inline

Collect distributed second-variational vectors into a global array.

Definition at line 345 of file k_point.hpp.

◆ gkvec()

template<typename T >
auto const & sirius::K_point< T >::gkvec ( ) const
inline

Definition at line 375 of file k_point.hpp.

◆ gkvec_sptr()

template<typename T >
auto sirius::K_point< T >::gkvec_sptr ( ) const
inline

Return shared pointer to gkvec object.

Definition at line 381 of file k_point.hpp.

◆ num_gkvec()

template<typename T >
int sirius::K_point< T >::num_gkvec ( ) const
inline

Total number of G+k vectors within the cutoff distance.

Definition at line 387 of file k_point.hpp.

◆ num_gkvec_loc()

template<typename T >
int sirius::K_point< T >::num_gkvec_loc ( ) const
inline

Local number of G+k vectors in case of flat distribution.

Definition at line 393 of file k_point.hpp.

◆ num_occupied_bands()

template<typename T >
int sirius::K_point< T >::num_occupied_bands ( int  ispn__ = -1) const
inline

Get the number of occupied bands for each spin channel.

Definition at line 399 of file k_point.hpp.

◆ band_energy() [1/2]

template<typename T >
double sirius::K_point< T >::band_energy ( int  j__,
int  ispn__ 
) const
inline

Get band energy.

Definition at line 413 of file k_point.hpp.

◆ band_energy() [2/2]

template<typename T >
void sirius::K_point< T >::band_energy ( int  j__,
int  ispn__,
double  e__ 
)
inline

Set band energy.

Definition at line 419 of file k_point.hpp.

◆ band_energies()

template<typename T >
auto sirius::K_point< T >::band_energies ( int  ispn__) const
inline

Definition at line 424 of file k_point.hpp.

◆ band_occupancy() [1/2]

template<typename T >
double sirius::K_point< T >::band_occupancy ( int  j__,
int  ispn__ 
) const
inline

Get band occupancy.

Definition at line 434 of file k_point.hpp.

◆ band_occupancy() [2/2]

template<typename T >
void sirius::K_point< T >::band_occupancy ( int  j__,
int  ispn__,
double  occ__ 
)
inline

Set band occupancy.

Definition at line 440 of file k_point.hpp.

◆ fv_eigen_value()

template<typename T >
double sirius::K_point< T >::fv_eigen_value ( int  i) const
inline

Definition at line 445 of file k_point.hpp.

◆ set_fv_eigen_values()

template<typename T >
void sirius::K_point< T >::set_fv_eigen_values ( double *  eval__)
inline

Definition at line 450 of file k_point.hpp.

◆ weight()

template<typename T >
double sirius::K_point< T >::weight ( ) const
inline

Return weight of k-point.

Definition at line 456 of file k_point.hpp.

◆ fv_states()

template<typename T >
auto & sirius::K_point< T >::fv_states ( )
inline

Definition at line 461 of file k_point.hpp.

◆ spinor_wave_functions() [1/2]

template<typename T >
wf::Wave_functions< T > const & sirius::K_point< T >::spinor_wave_functions ( ) const
inline

Definition at line 467 of file k_point.hpp.

◆ spinor_wave_functions() [2/2]

template<typename T >
wf::Wave_functions< T > & sirius::K_point< T >::spinor_wave_functions ( )
inline

Definition at line 473 of file k_point.hpp.

◆ spinor_wave_functions2()

template<typename T >
auto & sirius::K_point< T >::spinor_wave_functions2 ( )
inline

Definition at line 480 of file k_point.hpp.

◆ atomic_wave_functions_S() [1/2]

template<typename T >
auto const & sirius::K_point< T >::atomic_wave_functions_S ( ) const
inline

Return the initial atomic orbitals used to compute the hubbard wave functions. The S operator is applied on these functions.

Definition at line 489 of file k_point.hpp.

◆ atomic_wave_functions_S() [2/2]

template<typename T >
auto & sirius::K_point< T >::atomic_wave_functions_S ( )
inline

Definition at line 496 of file k_point.hpp.

◆ atomic_wave_functions() [1/2]

template<typename T >
auto const & sirius::K_point< T >::atomic_wave_functions ( ) const
inline

Return the initial atomic orbitals used to compute the hubbard wave functions.

Definition at line 502 of file k_point.hpp.

◆ atomic_wave_functions() [2/2]

template<typename T >
auto & sirius::K_point< T >::atomic_wave_functions ( )
inline

Return the initial atomic orbitals used to compute the hubbard wave functions.

Definition at line 509 of file k_point.hpp.

◆ hubbard_wave_functions_S()

template<typename T >
auto const & sirius::K_point< T >::hubbard_wave_functions_S ( ) const
inline

Return the actual hubbard wave functions used in the calculations.

The S operator is applied on these functions.

Definition at line 516 of file k_point.hpp.

◆ singular_components()

template<typename T >
auto & sirius::K_point< T >::singular_components ( )
inline

Definition at line 522 of file k_point.hpp.

◆ vk()

template<typename T >
auto sirius::K_point< T >::vk ( ) const
inline

Definition at line 527 of file k_point.hpp.

◆ gklo_basis_size()

template<typename T >
int sirius::K_point< T >::gklo_basis_size ( ) const
inline

Basis size of LAPW+lo method.

The total LAPW+lo basis size is equal to the sum of the number of LAPW functions and the total number of the local orbitals.

Definition at line 535 of file k_point.hpp.

◆ num_gkvec_row()

template<typename T >
int sirius::K_point< T >::num_gkvec_row ( ) const
inline

Local number of G+k vectors for each MPI rank in the row of the 2D MPI grid.

Definition at line 541 of file k_point.hpp.

◆ num_lo_row()

template<typename T >
int sirius::K_point< T >::num_lo_row ( ) const
inline

Local number of local orbitals for each MPI rank in the row of the 2D MPI grid.

Definition at line 547 of file k_point.hpp.

◆ gklo_basis_size_row()

template<typename T >
int sirius::K_point< T >::gklo_basis_size_row ( ) const
inline

Local number of basis functions for each MPI rank in the row of the 2D MPI grid.

Definition at line 553 of file k_point.hpp.

◆ num_gkvec_col()

template<typename T >
int sirius::K_point< T >::num_gkvec_col ( ) const
inline

Local number of G+k vectors for each MPI rank in the column of the 2D MPI grid.

Definition at line 559 of file k_point.hpp.

◆ num_lo_col()

template<typename T >
int sirius::K_point< T >::num_lo_col ( ) const
inline

Local number of local orbitals for each MPI rank in the column of the 2D MPI grid.

Definition at line 565 of file k_point.hpp.

◆ gklo_basis_size_col()

template<typename T >
int sirius::K_point< T >::gklo_basis_size_col ( ) const
inline

Local number of basis functions for each MPI rank in the column of the 2D MPI grid.

Definition at line 571 of file k_point.hpp.

◆ lo_basis_descriptor_col()

template<typename T >
lo_basis_descriptor const & sirius::K_point< T >::lo_basis_descriptor_col ( int  idx) const
inline

Definition at line 576 of file k_point.hpp.

◆ lo_basis_descriptor_row()

template<typename T >
lo_basis_descriptor const & sirius::K_point< T >::lo_basis_descriptor_row ( int  idx) const
inline

Definition at line 582 of file k_point.hpp.

◆ num_ranks_row()

template<typename T >
int sirius::K_point< T >::num_ranks_row ( ) const
inline

Definition at line 588 of file k_point.hpp.

◆ rank_row()

template<typename T >
int sirius::K_point< T >::rank_row ( ) const
inline

Definition at line 593 of file k_point.hpp.

◆ num_ranks_col()

template<typename T >
int sirius::K_point< T >::num_ranks_col ( ) const
inline

Definition at line 598 of file k_point.hpp.

◆ rank_col()

template<typename T >
int sirius::K_point< T >::rank_col ( ) const
inline

Definition at line 603 of file k_point.hpp.

◆ num_ranks()

template<typename T >
int sirius::K_point< T >::num_ranks ( ) const
inline

Number of MPI ranks for a given k-point.

Definition at line 609 of file k_point.hpp.

◆ rank()

template<typename T >
int sirius::K_point< T >::rank ( ) const
inline

Definition at line 614 of file k_point.hpp.

◆ num_atom_lo_cols()

template<typename T >
int sirius::K_point< T >::num_atom_lo_cols ( int  ia) const
inline

Return number of lo columns for a given atom.

Definition at line 620 of file k_point.hpp.

◆ lo_col()

template<typename T >
int sirius::K_point< T >::lo_col ( int  ia,
int  i 
) const
inline

Return local index (for the current MPI rank) of a column for a given atom and column index within an atom.

Definition at line 626 of file k_point.hpp.

◆ num_atom_lo_rows()

template<typename T >
int sirius::K_point< T >::num_atom_lo_rows ( int  ia) const
inline

Return number of lo rows for a given atom.

Definition at line 632 of file k_point.hpp.

◆ lo_row()

template<typename T >
int sirius::K_point< T >::lo_row ( int  ia,
int  i 
) const
inline

Return local index (for the current MPI rank) of a row for a given atom and row index within an atom.

Definition at line 638 of file k_point.hpp.

◆ fv_eigen_vectors()

template<typename T >
auto & sirius::K_point< T >::fv_eigen_vectors ( )
inline

Definition at line 643 of file k_point.hpp.

◆ fv_eigen_vectors_slab()

template<typename T >
auto & sirius::K_point< T >::fv_eigen_vectors_slab ( )
inline

Definition at line 648 of file k_point.hpp.

◆ sv_eigen_vectors()

template<typename T >
auto & sirius::K_point< T >::sv_eigen_vectors ( int  ispn)
inline

Definition at line 653 of file k_point.hpp.

◆ fd_eigen_vectors()

template<typename T >
auto & sirius::K_point< T >::fd_eigen_vectors ( )
inline

Definition at line 658 of file k_point.hpp.

◆ bypass_sv()

template<typename T >
void sirius::K_point< T >::bypass_sv ( )
inline

Definition at line 663 of file k_point.hpp.

◆ alm_coeffs_row()

template<typename T >
auto const & sirius::K_point< T >::alm_coeffs_row ( ) const
inline

Definition at line 668 of file k_point.hpp.

◆ alm_coeffs_col()

template<typename T >
auto const & sirius::K_point< T >::alm_coeffs_col ( ) const
inline

Definition at line 673 of file k_point.hpp.

◆ alm_coeffs_loc()

template<typename T >
auto const & sirius::K_point< T >::alm_coeffs_loc ( ) const
inline

Definition at line 678 of file k_point.hpp.

◆ comm()

template<typename T >
auto const & sirius::K_point< T >::comm ( ) const
inline

Definition at line 683 of file k_point.hpp.

◆ comm_row()

template<typename T >
auto const & sirius::K_point< T >::comm_row ( ) const
inline

Definition at line 688 of file k_point.hpp.

◆ comm_col()

template<typename T >
auto const & sirius::K_point< T >::comm_col ( ) const
inline

Definition at line 693 of file k_point.hpp.

◆ beta_projectors() [1/2]

template<typename T >
auto sirius::K_point< T >::beta_projectors ( ) -> Beta_projectors<T>&
inline

Definition at line 698 of file k_point.hpp.

◆ beta_projectors() [2/2]

template<typename T >
auto sirius::K_point< T >::beta_projectors ( ) const -> const Beta_projectors<T>&
inline

Definition at line 704 of file k_point.hpp.

◆ beta_projectors_row()

template<typename T >
auto & sirius::K_point< T >::beta_projectors_row ( )
inline

Definition at line 710 of file k_point.hpp.

◆ beta_projectors_col()

template<typename T >
auto & sirius::K_point< T >::beta_projectors_col ( )
inline

Definition at line 716 of file k_point.hpp.

◆ ctx()

template<typename T >
auto const & sirius::K_point< T >::ctx ( ) const
inline

Definition at line 722 of file k_point.hpp.

◆ out()

template<typename T >
std::ostream & sirius::K_point< T >::out ( int  level__) const
inline

Return stdout stream for high verbosity output.

This output will not be passed to a ctx_.out() stream.

Definition at line 729 of file k_point.hpp.

◆ spfft_transform()

template<typename T >
auto & sirius::K_point< T >::spfft_transform ( ) const
inline

Definition at line 738 of file k_point.hpp.

◆ gkvec_fft()

template<typename T >
auto const & sirius::K_point< T >::gkvec_fft ( ) const
inline

Definition at line 743 of file k_point.hpp.

◆ gkvec_fft_sptr()

template<typename T >
auto sirius::K_point< T >::gkvec_fft_sptr ( ) const
inline

Definition at line 748 of file k_point.hpp.

◆ gkvec_col()

template<typename T >
auto const & sirius::K_point< T >::gkvec_col ( ) const
inline

Definition at line 753 of file k_point.hpp.

◆ gkvec_row()

template<typename T >
auto const & sirius::K_point< T >::gkvec_row ( ) const
inline

Definition at line 758 of file k_point.hpp.

Friends And Related Function Documentation

◆ K_point_set

template<typename T >
friend class K_point_set
friend

Definition at line 199 of file k_point.hpp.

Member Data Documentation

◆ ctx_

template<typename T >
Simulation_context& sirius::K_point< T >::ctx_
private

Simulation context.

Definition at line 47 of file k_point.hpp.

◆ unit_cell_

template<typename T >
Unit_cell const& sirius::K_point< T >::unit_cell_
private

Unit cell object.

Definition at line 50 of file k_point.hpp.

◆ vk_

template<typename T >
r3::vector<double> sirius::K_point< T >::vk_
private

Fractional k-point coordinates.

Definition at line 53 of file k_point.hpp.

◆ weight_

template<typename T >
double sirius::K_point< T >::weight_ {1.0}
private

Weight of k-point.

Definition at line 56 of file k_point.hpp.

◆ comm_

template<typename T >
mpi::Communicator const& sirius::K_point< T >::comm_
private

Communicator for parallelization inside k-point.

This communicator is used to split G+k vectors and wave-functions.

Definition at line 60 of file k_point.hpp.

◆ gkvec_

template<typename T >
std::shared_ptr<fft::Gvec> sirius::K_point< T >::gkvec_
private

List of G-vectors with |G+k| < cutoff.

Definition at line 63 of file k_point.hpp.

◆ gkvec_row_

template<typename T >
std::shared_ptr<fft::Gvec> sirius::K_point< T >::gkvec_row_
private

Definition at line 65 of file k_point.hpp.

◆ gkvec_col_

template<typename T >
std::shared_ptr<fft::Gvec> sirius::K_point< T >::gkvec_col_
private

Definition at line 67 of file k_point.hpp.

◆ gkvec_partition_

template<typename T >
std::shared_ptr<fft::Gvec_fft> sirius::K_point< T >::gkvec_partition_
private

G-vector distribution for the FFT transformation.

Definition at line 70 of file k_point.hpp.

◆ spfft_transform_

template<typename T >
std::unique_ptr<fft::spfft_transform_type<T> > sirius::K_point< T >::spfft_transform_
mutableprivate

Definition at line 72 of file k_point.hpp.

◆ fv_eigen_values_

template<typename T >
sddk::mdarray<double, 1> sirius::K_point< T >::fv_eigen_values_
private

First-variational eigen values.

Definition at line 75 of file k_point.hpp.

◆ fv_eigen_vectors_

template<typename T >
la::dmatrix<std::complex<T> > sirius::K_point< T >::fv_eigen_vectors_
private

First-variational eigen vectors, distributed over 2D BLACS grid.

Definition at line 78 of file k_point.hpp.

◆ fv_eigen_vectors_slab_

template<typename T >
std::unique_ptr<wf::Wave_functions<T> > sirius::K_point< T >::fv_eigen_vectors_slab_
private

First-variational eigen vectors, distributed in slabs.

Definition at line 81 of file k_point.hpp.

◆ singular_components_

template<typename T >
std::unique_ptr<wf::Wave_functions<T> > sirius::K_point< T >::singular_components_
private

Lowest eigen-vectors of the LAPW overlap matrix with small aigen-values.

Definition at line 84 of file k_point.hpp.

◆ sv_eigen_vectors_

template<typename T >
std::array<la::dmatrix<std::complex<T> >, 2> sirius::K_point< T >::sv_eigen_vectors_
private

Second-variational eigen vectors.

Second-variational eigen-vectors are stored as one or two \( N_{fv} \times N_{fv} \) matrices in case of non-magnetic or collinear magnetic case or as a single \( 2 N_{fv} \times 2 N_{fv} \) matrix in case of general non-collinear magnetism.

Definition at line 90 of file k_point.hpp.

◆ fd_eigen_vectors_

template<typename T >
sddk::mdarray<std::complex<T>, 2> sirius::K_point< T >::fd_eigen_vectors_
private

Full-diagonalization eigen vectors.

Definition at line 93 of file k_point.hpp.

◆ fv_states_

template<typename T >
std::unique_ptr<wf::Wave_functions<T> > sirius::K_point< T >::fv_states_ {nullptr}
private

First-variational states.

Definition at line 96 of file k_point.hpp.

◆ spinor_wave_functions_

template<typename T >
std::unique_ptr<wf::Wave_functions<T> > sirius::K_point< T >::spinor_wave_functions_ {nullptr}
private

Two-component (spinor) wave functions describing the bands.

Definition at line 99 of file k_point.hpp.

◆ atomic_wave_functions_

template<typename T >
std::unique_ptr<wf::Wave_functions<T> > sirius::K_point< T >::atomic_wave_functions_ {nullptr}
private

Pseudopotential atmoic wave-functions (not orthogonalized).

Definition at line 102 of file k_point.hpp.

◆ atomic_wave_functions_S_

template<typename T >
std::unique_ptr<wf::Wave_functions<T> > sirius::K_point< T >::atomic_wave_functions_S_ {nullptr}
private

Pseudopotential atmoic wave-functions (not orthogonalized) with S-operator applied.

Definition at line 105 of file k_point.hpp.

◆ hubbard_wave_functions_

template<typename T >
std::unique_ptr<wf::Wave_functions<T> > sirius::K_point< T >::hubbard_wave_functions_ {nullptr}
private

Hubbard wave functions.

Definition at line 108 of file k_point.hpp.

◆ hubbard_wave_functions_S_

template<typename T >
std::unique_ptr<wf::Wave_functions<T> > sirius::K_point< T >::hubbard_wave_functions_S_ {nullptr}
private

Hubbard wave functions with S-operator applied.

Definition at line 111 of file k_point.hpp.

◆ band_occupancies_

template<typename T >
sddk::mdarray<double, 2> sirius::K_point< T >::band_occupancies_
private

Band occupation numbers.

Definition at line 114 of file k_point.hpp.

◆ band_energies_

template<typename T >
sddk::mdarray<double, 2> sirius::K_point< T >::band_energies_
private

Band energies.

Definition at line 117 of file k_point.hpp.

◆ alm_coeffs_row_

template<typename T >
std::unique_ptr<Matching_coefficients> sirius::K_point< T >::alm_coeffs_row_ {nullptr}
private

LAPW matching coefficients for the row G+k vectors.

Used to setup the distributed LAPW Hamiltonian and overlap matrices.

Definition at line 121 of file k_point.hpp.

◆ alm_coeffs_col_

template<typename T >
std::unique_ptr<Matching_coefficients> sirius::K_point< T >::alm_coeffs_col_ {nullptr}
private

LAPW matching coefficients for the column G+k vectors.

Used to setup the distributed LAPW Hamiltonian and overlap matrices.

Definition at line 125 of file k_point.hpp.

◆ alm_coeffs_loc_

template<typename T >
std::unique_ptr<Matching_coefficients> sirius::K_point< T >::alm_coeffs_loc_ {nullptr}
private

LAPW matching coefficients for the local set G+k vectors.

Definition at line 128 of file k_point.hpp.

◆ num_gkvec_row_

template<typename T >
int sirius::K_point< T >::num_gkvec_row_ {0}
private

Number of G+k vectors distributed along rows of MPI grid.

Definition at line 131 of file k_point.hpp.

◆ num_gkvec_col_

template<typename T >
int sirius::K_point< T >::num_gkvec_col_ {0}
private

Number of G+k vectors distributed along columns of MPI grid.

Definition at line 134 of file k_point.hpp.

◆ lo_basis_descriptors_row_

template<typename T >
std::vector<lo_basis_descriptor> sirius::K_point< T >::lo_basis_descriptors_row_
private

Basis descriptors distributed between rows of the 2D MPI grid.

This is a local array. Only MPI ranks belonging to the same column have identical copies of this array.

Definition at line 138 of file k_point.hpp.

◆ lo_basis_descriptors_col_

template<typename T >
std::vector<lo_basis_descriptor> sirius::K_point< T >::lo_basis_descriptors_col_
private

Basis descriptors distributed between columns of the 2D MPI grid.

This is a local array. Only MPI ranks belonging to the same row have identical copies of this array.

Definition at line 142 of file k_point.hpp.

◆ atom_lo_cols_

template<typename T >
std::vector<std::vector<int> > sirius::K_point< T >::atom_lo_cols_
private

List of columns of the Hamiltonian and overlap matrix lo block (local index) for a given atom.

Definition at line 145 of file k_point.hpp.

◆ atom_lo_rows_

template<typename T >
std::vector<std::vector<int> > sirius::K_point< T >::atom_lo_rows_
private

list of rows of the Hamiltonian and overlap matrix lo block (local index) for a given atom

Definition at line 148 of file k_point.hpp.

◆ zil_

template<typename T >
std::vector<std::complex<double> > sirius::K_point< T >::zil_
private

Imaginary unit to the power of l.

Definition at line 151 of file k_point.hpp.

◆ l_by_lm_

template<typename T >
std::vector<int> sirius::K_point< T >::l_by_lm_
private

Mapping between lm and l.

Definition at line 154 of file k_point.hpp.

◆ rank_col_

template<typename T >
int sirius::K_point< T >::rank_col_
private

Column rank of the processors of ScaLAPACK/ELPA diagonalization grid.

Definition at line 157 of file k_point.hpp.

◆ num_ranks_col_

template<typename T >
int sirius::K_point< T >::num_ranks_col_
private

Number of processors along the columns of the diagonalization grid.

Definition at line 160 of file k_point.hpp.

◆ rank_row_

template<typename T >
int sirius::K_point< T >::rank_row_
private

Row rank of the processors of ScaLAPACK/ELPA diagonalization grid.

Definition at line 163 of file k_point.hpp.

◆ num_ranks_row_

template<typename T >
int sirius::K_point< T >::num_ranks_row_
private

Number of processors along the rows of the diagonalization grid.

Definition at line 166 of file k_point.hpp.

◆ beta_projectors_

template<typename T >
std::unique_ptr<Beta_projectors<T> > sirius::K_point< T >::beta_projectors_ {nullptr}
private

Beta projectors for a local set of G+k vectors.

Definition at line 169 of file k_point.hpp.

◆ beta_projectors_row_

template<typename T >
std::unique_ptr<Beta_projectors<T> > sirius::K_point< T >::beta_projectors_row_ {nullptr}
private

Beta projectors for row G+k vectors.

Used to setup the full Hamiltonian in PP-PW case (for verification purpose only)

Definition at line 173 of file k_point.hpp.

◆ beta_projectors_col_

template<typename T >
std::unique_ptr<Beta_projectors<T> > sirius::K_point< T >::beta_projectors_col_ {nullptr}
private

Beta projectors for column G+k vectors.

Used to setup the full Hamiltonian in PP-PW case (for verification purpose only)

Definition at line 177 of file k_point.hpp.

◆ comm_row_

template<typename T >
mpi::Communicator const& sirius::K_point< T >::comm_row_
private

Communicator between(!!) rows.

Definition at line 180 of file k_point.hpp.

◆ comm_col_

template<typename T >
mpi::Communicator const& sirius::K_point< T >::comm_col_
private

Communicator between(!!) columns.

Definition at line 183 of file k_point.hpp.

◆ ispn_map_

template<typename T >
std::array<int, 2> sirius::K_point< T >::ispn_map_ {0, -1}
private

Definition at line 185 of file k_point.hpp.


The documentation for this class was generated from the following files: