SIRIUS 7.5.0
Electronic structure library and applications
Public Member Functions | Private Member Functions | Private Attributes | List of all members
sirius::Density Class Reference

Generate charge density and magnetization from occupied spinor wave-functions. More...

#include <density.hpp>

Inherits sirius::Field4D.

Public Member Functions

 Density (Simulation_context &ctx__)
 Constructor. More...
 
void update ()
 Update internal parameters after a change of lattice vectors or atomic coordinates. More...
 
double core_leakage () const
 Find the total leakage of the core states out of the muffin-tins. More...
 
void initial_density ()
 Generate initial charge density and magnetization. More...
 
void initial_density_pseudo ()
 
void initial_density_full_pot ()
 
void normalize ()
 
bool check_num_electrons () const
 Check total density for the correct number of electrons. More...
 
template<typename T >
void generate (K_point_set const &ks__, bool symmetrize__, bool add_core__, bool transform_to_rg__)
 Generate full charge density (valence + core) and magnetization from the wave functions. More...
 
template<typename T >
void generate_valence (K_point_set const &ks__)
 Generate valence charge density and magnetization from the wave functions. More...
 
void augment ()
 Add augmentation charge Q(r). More...
 
sddk::mdarray< std::complex< double >, 2 > generate_rho_aug () const
 Generate augmentation charge density. More...
 
double core_leakage (int ic) const
 Return core leakage for a specific atom symmetry class. More...
 
auto const & rho () const
 Return const reference to charge density (scalar functions). More...
 
auto & rho ()
 Return charge density (scalar functions). More...
 
auto & mag (int i)
 
auto const & mag (int i) const
 
auto & rho_pseudo_core ()
 
auto const & rho_pseudo_core () const
 
auto const & density_mt (atom_index_t::local ialoc__) const
 
void generate_paw_density ()
 Generate \( n_1 \) and \( \tilde{n}_1 \) in lm components. More...
 
auto paw_ae_density (int ia__) const
 Return list of pointers to all-electron PAW density function for a given local index of atom with PAW potential. More...
 
auto paw_ps_density (int ia__) const
 Return list of pointers to pseudo PAW density function for a given local index of atom with PAW potential. More...
 
void mixer_input ()
 
void mixer_output ()
 
void mixer_init (config_t::mixer_t const &mixer_cfg__)
 Initialize density mixer. More...
 
double mix ()
 Mix new density. More...
 
auto const & density_matrix (int ia__) const
 
auto & density_matrix (int ia__)
 
sddk::mdarray< double, 3 > density_matrix_aux (Atom_type const &atom_type__) const
 Return density matrix for all atoms of a given type in auxiliary form. More...
 
sddk::mdarray< double, 2 > density_matrix_aux (typename atom_index_t::global ia__) const
 
sddk::mdarray< double, 2 > compute_atomic_mag_mom () const
 Calculate approximate atomic magnetic moments in case of PP-PW. More...
 
std::tuple< std::array< double, 3 >, std::array< double, 3 >, std::vector< std::array< double, 3 > > > get_magnetisation () const
 Get total magnetization and also contributions from interstitial and muffin-tin parts. More...
 
void print_info (std::ostream &out__) const
 
Occupation_matrix const & occupation_matrix () const
 
Occupation_matrixoccupation_matrix ()
 
auto const & paw_density () const
 
void check_density_continuity_at_mt ()
 Check density at MT boundary. More...
 
void save (std::string name__) const
 
void load (std::string name__)
 
void save_to_xsf ()
 
void save_to_ted ()
 
void save_to_xdmf ()
 
- Public Member Functions inherited from sirius::Field4D
 Field4D (Simulation_context &ctx__, lmax_t lmax__, std::array< periodic_function_ptr_t< double > const *, 4 > ptr__={nullptr, nullptr, nullptr, nullptr})
 Constructor. More...
 
auto & scalar ()
 Return scalar part of the field. More...
 
auto const & scalar () const
 Return scalar part of the field. More...
 
auto & vector (int i)
 Return component of the vector part of the field. More...
 
auto const & vector (int i) const
 Return component of the vector part of the field. More...
 
auto & component (int i)
 
auto & component_raise (int i)
 Throws error in case of invalid access. More...
 
auto const & component (int i) const
 
void zero ()
 
void fft_transform (int direction__)
 
auto & ctx ()
 
auto const & ctx () const
 
auto mt_components ()
 
auto pw_components ()
 

Private Member Functions

void generate_paw_density (paw_atom_index_t::local iapaw__)
 Generate atomic densities in the case of PAW. More...
 
void init_density_matrix_for_paw ()
 Initialize PAW density matrix. More...
 
template<int num_mag_dims>
void reduce_density_matrix (Atom_type const &atom_type__, sddk::mdarray< std::complex< double >, 3 > const &zdens__, sddk::mdarray< double, 3 > &mt_density_matrix__)
 Reduce complex density matrix over magnetic quantum numbers. More...
 
template<typename T , typename F >
void add_k_point_contribution_dm (K_point< T > &kp__, sddk::mdarray< std::complex< double >, 4 > &density_matrix__)
 Add k-point contribution to the density matrix in the canonical form. More...
 
template<typename T >
void add_k_point_contribution_rg (K_point< T > *kp__, std::array< wf::Wave_functions_fft< T >, 2 > &wf_fft__)
 Add k-point contribution to the density and magnetization defined on the regular FFT grid. More...
 
void generate_valence_mt ()
 Generate valence density in the muffin-tins. More...
 
void generate_core_charge_density ()
 Generate charge density of core states. More...
 
void generate_pseudo_core_charge_density ()
 

Private Attributes

Unit_cellunit_cell_
 Alias to ctx_.unit_cell() More...
 
std::unique_ptr< density_matrix_tdensity_matrix_
 Density matrix for all atoms. More...
 
std::unique_ptr< PAW_density< double > > paw_density_
 Local fraction of atoms with PAW correction. More...
 
std::unique_ptr< Occupation_matrixoccupation_matrix_
 Occupation matrix of the LDA+U method. More...
 
std::array< std::unique_ptr< Smooth_periodic_function< double > >, 4 > rho_mag_coarse_
 Density and magnetization on the coarse FFT mesh. More...
 
std::unique_ptr< Smooth_periodic_function< double > > rho_pseudo_core_ {nullptr}
 Pointer to pseudo core charge density. More...
 
std::vector< int > l_by_lm_
 Fast mapping between composite lm index and corresponding orbital quantum number. More...
 
std::unique_ptr< mixer::Mixer< Periodic_function< double >, Periodic_function< double >, Periodic_function< double >, Periodic_function< double >, density_matrix_t, PAW_density< double >, Hubbard_matrix > > mixer_
 Density mixer. More...
 

Additional Inherited Members

- Protected Attributes inherited from sirius::Field4D
Simulation_contextctx_
 

Detailed Description

Generate charge density and magnetization from occupied spinor wave-functions.

Let's start from the definition of the complex density matrix:

\[ \rho_{\sigma' \sigma}({\bf r}) = \sum_{j{\bf k}} n_{j{\bf k}} \Psi_{j{\bf k}}^{\sigma*}({\bf r}) \Psi_{j{\bf k}}^{\sigma'}({\bf r}) = \frac{1}{2} \left( \begin{array}{cc} \rho({\bf r})+m_z({\bf r}) & m_x({\bf r})-im_y({\bf r}) \\ m_x({\bf r})+im_y({\bf r}) & \rho({\bf r})-m_z({\bf r}) \end{array} \right) \]

We notice that the diagonal components of the density matrix are actually real and the off-diagonal components are expressed trough two independent functions \( m_x({\bf r}) \) and \( m_y({\bf r}) \). Having this in mind we will work with a slightly different object, namely a real density matrix, defined as a 1-, 2- or 4-dimensional (depending on the number of magnetic components) vector with the following elements:

At this point it is straightforward to compute the density and magnetization in the interstitial (see Density::add_k_point_contribution_rg()). The muffin-tin part of the density and magnetization is obtained in a slighlty more complicated way. Recall the expansion of spinor wave-functions inside the muffin-tin \( \alpha \):

\[ \Psi_{j{\bf k}}^{\sigma}({\bf r}) = \sum_{\xi}^{N_{\xi}^{\alpha}} {S_{\xi}^{\sigma j {\bf k},\alpha}} f_{\ell_{\xi} \lambda_{\xi}}^{\alpha}(r)Y_{\ell_{\xi}m_{\xi}}(\hat {\bf r}) \]

which we insert into expression for the complex density matrix:

\[ \rho_{\sigma' \sigma}({\bf r}) = \sum_{j{\bf k}} n_{j{\bf k}} \sum_{\xi}^{N_{\xi}^{\alpha}} S_{\xi}^{\sigma j {\bf k},\alpha*} f_{\ell_{\xi} \lambda_{\xi}}^{\alpha}(r) Y_{\ell_{\xi}m_{\xi}}^{*}(\hat {\bf r}) \sum_{\xi'}^{N_{\xi'}^{\alpha}} S_{\xi'}^{\sigma' j{\bf k},\alpha} f_{\ell_{\xi'} \lambda_{\xi'}}^{\alpha}(r)Y_{\ell_{\xi'}m_{\xi'}}(\hat {\bf r}) \]

First, we eliminate a sum over bands and k-points by forming an auxiliary density tensor:

\[ D_{\xi \sigma, \xi' \sigma'}^{\alpha} = \sum_{j{\bf k}} n_{j{\bf k}} S_{\xi}^{\sigma j {\bf k},\alpha*} S_{\xi'}^{\sigma' j {\bf k},\alpha} \]

The expression for complex density matrix simplifies to:

\[ \rho_{\sigma' \sigma}({\bf r}) = \sum_{\xi \xi'} D_{\xi \sigma, \xi' \sigma'}^{\alpha} f_{\ell_{\xi} \lambda_{\xi}}^{\alpha}(r)Y_{\ell_{\xi}m_{\xi}}^{*}(\hat {\bf r}) f_{\ell_{\xi'} \lambda_{\xi'}}^{\alpha}(r)Y_{\ell_{\xi'}m_{\xi'}}(\hat {\bf r}) \]

Now we can switch to the real density matrix and write its' expansion in real spherical harmonics. Let's take non-magnetic case as an example:

\[ \rho({\bf r}) = \sum_{\xi \xi'} D_{\xi \xi'}^{\alpha} f_{\ell_{\xi} \lambda_{\xi}}^{\alpha}(r)Y_{\ell_{\xi}m_{\xi}}^{*}(\hat {\bf r}) f_{\ell_{\xi'} \lambda_{\xi'}}^{\alpha}(r)Y_{\ell_{\xi'}m_{\xi'}}(\hat {\bf r}) = \sum_{\ell_3 m_3} \rho_{\ell_3 m_3}^{\alpha}(r) R_{\ell_3 m_3}(\hat {\bf r}) \]

where

\[ \rho_{\ell_3 m_3}^{\alpha}(r) = \sum_{\xi \xi'} D_{\xi \xi'}^{\alpha} f_{\ell_{\xi} \lambda_{\xi}}^{\alpha}(r) f_{\ell_{\xi'} \lambda_{\xi'}}^{\alpha}(r) \langle Y_{\ell_{\xi}m_{\xi}} | R_{\ell_3 m_3} | Y_{\ell_{\xi'}m_{\xi'}} \rangle \]

We are almost done. Now it is time to switch to the full index notation \( \xi \rightarrow \{ \ell \lambda m \} \) and sum over m and m' indices:

\[ \rho_{\ell_3 m_3}^{\alpha}(r) = \sum_{\ell \lambda, \ell' \lambda'} f_{\ell \lambda}^{\alpha}(r) f_{\ell' \lambda'}^{\alpha}(r) d_{\ell \lambda, \ell' \lambda', \ell_3 m_3}^{\alpha} \]

where

\[ d_{\ell \lambda, \ell' \lambda', \ell_3 m_3}^{\alpha} = \sum_{mm'} D_{\ell \lambda m, \ell' \lambda' m'}^{\alpha} \langle Y_{\ell m} | R_{\ell_3 m_3} | Y_{\ell' m'} \rangle \]

This is our final answer: radial components of density and magnetization are expressed as a linear combination of quadratic forms in radial functions.

Note
density and potential are allocated as global function because it's easier to load and save them.
In case of full-potential calculation valence + core electron charge density is computed.
In tcase of pseudopotential valence charge density is computed.

Definition at line 213 of file density.hpp.

Constructor & Destructor Documentation

◆ Density()

sirius::Density::Density ( Simulation_context ctx__)

Constructor.

Definition at line 81 of file density.cpp.

Member Function Documentation

◆ generate_paw_density() [1/2]

void sirius::Density::generate_paw_density ( paw_atom_index_t::local  iapaw__)
private

Generate atomic densities in the case of PAW.

Definition at line 482 of file density.cpp.

◆ init_density_matrix_for_paw()

void sirius::Density::init_density_matrix_for_paw ( )
private

Initialize PAW density matrix.

Definition at line 436 of file density.cpp.

◆ reduce_density_matrix()

template<int num_mag_dims>
void sirius::Density::reduce_density_matrix ( Atom_type const &  atom_type__,
sddk::mdarray< std::complex< double >, 3 > const &  zdens__,
sddk::mdarray< double, 3 > &  mt_density_matrix__ 
)
private

Reduce complex density matrix over magnetic quantum numbers.

The following operation is performed:

\[ n_{\ell \lambda, \ell' \lambda', \ell_3 m_3}^{\alpha} = \sum_{mm'} D_{\ell \lambda m, \ell' \lambda' m'}^{\alpha} \langle Y_{\ell m} | R_{\ell_3 m_3} | Y_{\ell' m'} \rangle \]

Definition at line 1541 of file density.cpp.

◆ add_k_point_contribution_dm()

template<typename T , typename F >
void sirius::Density::add_k_point_contribution_dm ( K_point< T > &  kp__,
sddk::mdarray< std::complex< double >, 4 > &  density_matrix__ 
)
private

Add k-point contribution to the density matrix in the canonical form.

In case of full-potential LAPW complex density matrix has the following expression:

\[ d_{\xi \sigma, \xi' \sigma'}^{\alpha} = \sum_{j{\bf k}} n_{j{\bf k}} S_{\xi}^{\sigma j {\bf k},\alpha*} S_{\xi'}^{\sigma' j {\bf k},\alpha} \]

where \( S_{\xi}^{\sigma j {\bf k},\alpha} \) are the expansion coefficients of spinor wave functions inside muffin-tin spheres.

In case of LDA+U the occupation matrix is also computed. It has the following expression:

\[ n_{\ell,mm'}^{\sigma \sigma'} = \sum_{i {\bf k}}^{occ} \int_{0}^{R_{MT}} r^2 dr \Psi_{\ell m}^{i{\bf k}\sigma *}({\bf r}) \Psi_{\ell m'}^{i{\bf k}\sigma'}({\bf r}) \]

In case of ultrasoft pseudopotential the following density matrix has to be computed for each atom:

\[ d_{\xi \xi'}^{\alpha} = \langle \beta_{\xi}^{\alpha} | \hat N | \beta_{\xi'}^{\alpha} \rangle = \sum_{j {\bf k}} \langle \beta_{\xi}^{\alpha} | \Psi_{j{\bf k}} \rangle n_{j{\bf k}} \langle \Psi_{j{\bf k}} | \beta_{\xi'}^{\alpha} \rangle \]

Here \( \hat N = \sum_{j{\bf k}} | \Psi_{j{\bf k}} \rangle n_{j{\bf k}} \langle \Psi_{j{\bf k}} | \) is the occupancy operator written in spectral representation.

Template Parameters
TPrecision type of wave-functions.
FType of the wave-functions inner product (used in pp-pw).

◆ add_k_point_contribution_rg()

template<typename T >
void sirius::Density::add_k_point_contribution_rg ( K_point< T > *  kp__,
std::array< wf::Wave_functions_fft< T >, 2 > &  wf_fft__ 
)
private

Add k-point contribution to the density and magnetization defined on the regular FFT grid.

Definition at line 676 of file density.cpp.

◆ generate_valence_mt()

void sirius::Density::generate_valence_mt ( )
private

Generate valence density in the muffin-tins.

Definition at line 1580 of file density.cpp.

◆ generate_core_charge_density()

void sirius::Density::generate_core_charge_density ( )
inlineprivate

Generate charge density of core states.

Definition at line 311 of file density.hpp.

◆ generate_pseudo_core_charge_density()

void sirius::Density::generate_pseudo_core_charge_density ( )
inlineprivate

Definition at line 327 of file density.hpp.

◆ update()

void sirius::Density::update ( )

Update internal parameters after a change of lattice vectors or atomic coordinates.

Definition at line 120 of file density.cpp.

◆ core_leakage() [1/2]

double sirius::Density::core_leakage ( ) const

Find the total leakage of the core states out of the muffin-tins.

Definition at line 138 of file density.cpp.

◆ initial_density()

void sirius::Density::initial_density ( )

Generate initial charge density and magnetization.

Definition at line 148 of file density.cpp.

◆ initial_density_pseudo()

void sirius::Density::initial_density_pseudo ( )

Definition at line 173 of file density.cpp.

◆ initial_density_full_pot()

void sirius::Density::initial_density_full_pot ( )

Definition at line 270 of file density.cpp.

◆ normalize()

void sirius::Density::normalize ( )

Definition at line 1035 of file density.cpp.

◆ check_num_electrons()

bool sirius::Density::check_num_electrons ( ) const

Check total density for the correct number of electrons.

Definition at line 1057 of file density.cpp.

◆ generate()

template<typename T >
template void sirius::Density::generate< double > ( K_point_set const &  ks__,
bool  symmetrize__,
bool  add_core__,
bool  transform_to_rg__ 
)

Generate full charge density (valence + core) and magnetization from the wave functions.

This function calls generate_valence() and then in case of full-potential LAPW method adds a core density to get the full charge density of the system. Density is generated in spectral representation, i.e. plane-wave coefficients in the interstitial and spherical harmonic components in the muffin-tins.

Definition at line 1088 of file density.cpp.

◆ generate_valence()

template<typename T >
void sirius::Density::generate_valence ( K_point_set const &  ks__)

Generate valence charge density and magnetization from the wave functions.

The interstitial density is generated on the coarse FFT grid and then transformed to the PW domain. After symmetrization and mixing and before the generation of the XC potential density is transformted to the real-space domain and checked for the number of electrons.

Definition at line 1230 of file density.cpp.

◆ augment()

void sirius::Density::augment ( )

Add augmentation charge Q(r).

Restore valence density by adding the Q-operator constribution. The following term is added to the valence density, generated by the pseudo wave-functions:

\[ \tilde \rho({\bf G}) = \sum_{\alpha} \sum_{\xi \xi'} d_{\xi \xi'}^{\alpha} Q_{\xi' \xi}^{\alpha}({\bf G}) \]

Plane-wave coefficients of the Q-operator for a given atom \( \alpha \) can be obtained from the corresponding coefficients of the Q-operator for a given atom type A:

\[ Q_{\xi' \xi}^{\alpha(A)}({\bf G}) = e^{-i{\bf G}\tau_{\alpha(A)}} Q_{\xi' \xi}^{A}({\bf G}) \]

We use this property to split the sum over atoms into sum over atom types and inner sum over atoms of the same type:

\[ \tilde \rho({\bf G}) = \sum_{A} \sum_{\xi \xi'} Q_{\xi' \xi}^{A}({\bf G}) \sum_{\alpha(A)} d_{\xi \xi'}^{\alpha(A)} e^{-i{\bf G}\tau_{\alpha(A)}} = \sum_{A} \sum_{\xi \xi'} Q_{\xi' \xi}^{A}({\bf G}) d_{\xi \xi'}^{A}({\bf G}) \]

where

\[ d_{\xi \xi'}^{A}({\bf G}) = \sum_{\alpha(A)} d_{\xi \xi'}^{\alpha(A)} e^{-i{\bf G}\tau_{\alpha(A)}} \]

Definition at line 1209 of file density.cpp.

◆ generate_rho_aug()

sddk::mdarray< std::complex< double >, 2 > sirius::Density::generate_rho_aug ( ) const

Generate augmentation charge density.

Definition at line 1374 of file density.cpp.

◆ core_leakage() [2/2]

double sirius::Density::core_leakage ( int  ic) const
inline

Return core leakage for a specific atom symmetry class.

Definition at line 410 of file density.hpp.

◆ rho() [1/2]

auto const & sirius::Density::rho ( ) const
inline

Return const reference to charge density (scalar functions).

Definition at line 416 of file density.hpp.

◆ rho() [2/2]

auto & sirius::Density::rho ( )
inline

Return charge density (scalar functions).

Definition at line 422 of file density.hpp.

◆ mag() [1/2]

auto & sirius::Density::mag ( int  i)
inline

Definition at line 427 of file density.hpp.

◆ mag() [2/2]

auto const & sirius::Density::mag ( int  i) const
inline

Definition at line 432 of file density.hpp.

◆ rho_pseudo_core() [1/2]

auto & sirius::Density::rho_pseudo_core ( )
inline

Definition at line 437 of file density.hpp.

◆ rho_pseudo_core() [2/2]

auto const & sirius::Density::rho_pseudo_core ( ) const
inline

Definition at line 442 of file density.hpp.

◆ density_mt()

auto const & sirius::Density::density_mt ( atom_index_t::local  ialoc__) const
inline

Definition at line 447 of file density.hpp.

◆ generate_paw_density() [2/2]

void sirius::Density::generate_paw_density ( )

Generate \( n_1 \) and \( \tilde{n}_1 \) in lm components.

Definition at line 552 of file density.cpp.

◆ paw_ae_density()

auto sirius::Density::paw_ae_density ( int  ia__) const
inline

Return list of pointers to all-electron PAW density function for a given local index of atom with PAW potential.

Definition at line 457 of file density.hpp.

◆ paw_ps_density()

auto sirius::Density::paw_ps_density ( int  ia__) const
inline

Return list of pointers to pseudo PAW density function for a given local index of atom with PAW potential.

Definition at line 467 of file density.hpp.

◆ mixer_input()

void sirius::Density::mixer_input ( )

Definition at line 1880 of file density.cpp.

◆ mixer_output()

void sirius::Density::mixer_output ( )

Definition at line 1905 of file density.cpp.

◆ mixer_init()

void sirius::Density::mixer_init ( config_t::mixer_t const &  mixer_cfg__)

Initialize density mixer.

Definition at line 1830 of file density.cpp.

◆ mix()

double sirius::Density::mix ( )

Mix new density.

Definition at line 1933 of file density.cpp.

◆ density_matrix() [1/2]

auto const & sirius::Density::density_matrix ( int  ia__) const
inline

Definition at line 486 of file density.hpp.

◆ density_matrix() [2/2]

auto & sirius::Density::density_matrix ( int  ia__)
inline

Definition at line 491 of file density.hpp.

◆ density_matrix_aux() [1/2]

sddk::mdarray< double, 3 > sirius::Density::density_matrix_aux ( Atom_type const &  atom_type__) const

Return density matrix for all atoms of a given type in auxiliary form.

Definition at line 1810 of file density.cpp.

◆ density_matrix_aux() [2/2]

sddk::mdarray< double, 2 > sirius::Density::density_matrix_aux ( typename atom_index_t::global  ia__) const

Definition at line 1780 of file density.cpp.

◆ compute_atomic_mag_mom()

sddk::mdarray< double, 2 > sirius::Density::compute_atomic_mag_mom ( ) const

Calculate approximate atomic magnetic moments in case of PP-PW.

Definition at line 1718 of file density.cpp.

◆ get_magnetisation()

std::tuple< std::array< double, 3 >, std::array< double, 3 >, std::vector< std::array< double, 3 > > > sirius::Density::get_magnetisation ( ) const

Get total magnetization and also contributions from interstitial and muffin-tin parts.

In case of PP-PW there are no real muffin-tins. Instead, a value of magnetization inside atomic sphere with some chosen radius is returned.

Definition at line 1746 of file density.cpp.

◆ print_info()

void sirius::Density::print_info ( std::ostream &  out__) const

Definition at line 1944 of file density.cpp.

◆ occupation_matrix() [1/2]

Occupation_matrix const & sirius::Density::occupation_matrix ( ) const
inline

Definition at line 516 of file density.hpp.

◆ occupation_matrix() [2/2]

Occupation_matrix & sirius::Density::occupation_matrix ( )
inline

Definition at line 521 of file density.hpp.

◆ paw_density()

auto const & sirius::Density::paw_density ( ) const
inline

Definition at line 526 of file density.hpp.

◆ check_density_continuity_at_mt()

void sirius::Density::check_density_continuity_at_mt ( )
inline

Check density at MT boundary.

Definition at line 532 of file density.hpp.

◆ save()

void sirius::Density::save ( std::string  name__) const
inline

Definition at line 568 of file density.hpp.

◆ load()

void sirius::Density::load ( std::string  name__)
inline

Definition at line 579 of file density.hpp.

◆ save_to_xsf()

void sirius::Density::save_to_xsf ( )
inline

Definition at line 599 of file density.hpp.

◆ save_to_ted()

void sirius::Density::save_to_ted ( )
inline

Definition at line 624 of file density.hpp.

◆ save_to_xdmf()

void sirius::Density::save_to_xdmf ( )
inline

Definition at line 654 of file density.hpp.

Member Data Documentation

◆ unit_cell_

Unit_cell& sirius::Density::unit_cell_
private

Alias to ctx_.unit_cell()

Definition at line 217 of file density.hpp.

◆ density_matrix_

std::unique_ptr<density_matrix_t> sirius::Density::density_matrix_
private

Density matrix for all atoms.

This is a global matrix, meaning that each MPI rank holds the full copy. This simplifies the symmetrization.

Definition at line 221 of file density.hpp.

◆ paw_density_

std::unique_ptr<PAW_density<double> > sirius::Density::paw_density_
private

Local fraction of atoms with PAW correction.

Definition at line 224 of file density.hpp.

◆ occupation_matrix_

std::unique_ptr<Occupation_matrix> sirius::Density::occupation_matrix_
private

Occupation matrix of the LDA+U method.

Definition at line 227 of file density.hpp.

◆ rho_mag_coarse_

std::array<std::unique_ptr<Smooth_periodic_function<double> >, 4> sirius::Density::rho_mag_coarse_
private

Density and magnetization on the coarse FFT mesh.

Coarse FFT grid is enough to generate density and magnetization from the wave-functions. The components of the rho_mag_coarse vector have the following order: \( \{\rho({\bf r}), m_z({\bf r}), m_x({\bf r}), m_y({\bf r}) \} \).

Definition at line 234 of file density.hpp.

◆ rho_pseudo_core_

std::unique_ptr<Smooth_periodic_function<double> > sirius::Density::rho_pseudo_core_ {nullptr}
private

Pointer to pseudo core charge density.

In the case of pseudopotential we need to know the non-linear core correction to the exchange-correlation energy which is introduced trough the pseudo core density: \( E_{xc}[\rho_{val} + \rho_{core}] \). The 'pseudo' reflects the fact that this density integrated does not reproduce the total number of core elctrons.

Definition at line 242 of file density.hpp.

◆ l_by_lm_

std::vector<int> sirius::Density::l_by_lm_
private

Fast mapping between composite lm index and corresponding orbital quantum number.

Definition at line 245 of file density.hpp.

◆ mixer_

std::unique_ptr<mixer::Mixer<Periodic_function<double>, Periodic_function<double>, Periodic_function<double>, Periodic_function<double>, density_matrix_t, PAW_density<double>, Hubbard_matrix> > sirius::Density::mixer_
private

Density mixer.

Mix the following objects: density, x-,y-,z-components of magnetisation, density matrix and PAW density of atoms.

Definition at line 252 of file density.hpp.


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