SIRIUS 7.5.0
Electronic structure library and applications
|
#include <stress.hpp>
Public Member Functions | |
Stress (Simulation_context &ctx__, Density &density__, Potential &potential__, K_point_set &kset__) | |
r3::matrix< double > | calc_stress_vloc () |
Local potential contribution to stress. More... | |
r3::matrix< double > | stress_vloc () const |
r3::matrix< double > | calc_stress_har () |
Hartree energy contribution to stress. More... | |
r3::matrix< double > | stress_har () const |
r3::matrix< double > | calc_stress_ewald () |
Ewald energy contribution to stress. More... | |
r3::matrix< double > | stress_ewald () const |
template<typename T > | |
void | calc_stress_kin_aux () |
Kinetic energy contribution to stress. More... | |
r3::matrix< double > | calc_stress_kin () |
r3::matrix< double > | stress_kin () const |
r3::matrix< double > | calc_stress_nonloc () |
r3::matrix< double > | stress_nonloc () const |
r3::matrix< double > | calc_stress_us () |
Contribution to the stress tensor from the augmentation operator. More... | |
auto | stress_us () const |
auto | stress_us_nl () const |
r3::matrix< double > | calc_stress_xc () |
XC contribution to stress. More... | |
auto | stress_xc () const |
r3::matrix< double > | calc_stress_core () |
Non-linear core correction to stress tensor. More... | |
auto | stress_core () const |
r3::matrix< double > | calc_stress_hubbard () |
auto | stress_hubbard () const |
r3::matrix< double > | calc_stress_total () |
auto | stress_total () const |
void | print_info (std::ostream &out__, int verbosity__) const |
Private Member Functions | |
template<typename T , typename F > | |
void | calc_stress_nonloc_aux () |
Non-local contribution to stress. More... | |
Private Attributes | |
Simulation_context & | ctx_ |
Density const & | density_ |
Potential & | potential_ |
K_point_set & | kset_ |
r3::matrix< double > | stress_kin_ |
r3::matrix< double > | stress_har_ |
r3::matrix< double > | stress_ewald_ |
r3::matrix< double > | stress_vloc_ |
r3::matrix< double > | stress_nonloc_ |
r3::matrix< double > | stress_us_ |
r3::matrix< double > | stress_xc_ |
r3::matrix< double > | stress_core_ |
r3::matrix< double > | stress_hubbard_ |
r3::matrix< double > | stress_total_ |
Stress tensor.
The following referenceces were particularly useful in the derivation of the stress tensor components:
Stress tensor describes a reaction of crystall to a strain:
\[ \sigma_{\mu \nu} = \frac{1}{\Omega} \frac{\partial E}{\partial \varepsilon_{\mu \nu}} \]
where \( \varepsilon_{\mu \nu} \) is a symmetric strain tensor which describes the infinitesimal deformation of a crystal:
\[ r_{\mu} \rightarrow \sum_{\nu} (\delta_{\mu \nu} + \varepsilon_{\mu \nu}) r_{\nu} \]
The following expressions are helpful:
\[ \frac{\partial r_{\tau}}{\partial \varepsilon_{\mu \nu}} = \delta_{\tau \mu} r_{\nu} \]
\[ \frac{\partial |{\bf r}|}{\partial \varepsilon_{\mu \nu}} = \sum_{\tau} \frac{\partial |{\bf r}|}{\partial r_{\tau}} \frac{\partial r_{\tau}}{\partial \varepsilon_{\mu \nu}} = \sum_{\tau} \frac{ r_{\tau}}{|{\bf r}|} \delta_{\tau \mu} r_{\nu} = \frac{r_{\mu} r_{\nu}}{|{\bf r}|} \]
\[ \frac{\partial \Omega}{\partial \varepsilon_{\mu \nu}} = \delta_{\mu \nu} \Omega \]
\[ \frac{\partial}{\partial \varepsilon_{\mu \nu}} \frac{1}{\sqrt{\Omega}} = -\frac{1}{2}\frac{1}{\sqrt{\Omega}} \delta_{\mu \nu} \]
\[ \frac{\partial G_{\tau}}{\partial \varepsilon_{\mu \nu}} = -\delta_{\tau \nu} G_{\mu} \]
\[ \frac{\partial |{\bf G}|}{\partial \varepsilon_{\mu \nu}} = \sum_{\tau} \frac{\partial |{\bf G}|}{\partial G_{\tau}} \frac{\partial G_{\tau}}{\partial \varepsilon_{\mu \nu}} = -\frac{1}{|{\bf G}|}G_{\nu}G_{\mu} \]
In the derivation of the stress tensor contributions it is important to know which variables are invariant under lattice distortion. This are:\[ \psi({\bf G}) = \frac{1}{\sqrt{\Omega}} \int e^{-i {\bf G}{\bf r}} \psi({\bf r}) d{\bf r} \]
\[ \tilde \rho({\bf G}) = \int e^{-i {\bf G}{\bf r}} \sum_{j} |\psi_j({\bf r})|^{2} d{\bf r} \]
Definition at line 96 of file stress.hpp.
|
inline |
Definition at line 222 of file stress.hpp.
|
private |
Non-local contribution to stress.
Energy contribution from the non-local part of pseudopotential:
\[ E^{nl} = \sum_{i{\bf k}} \sum_{\alpha}\sum_{\xi \xi'}P_{\xi}^{\alpha,i{\bf k}} D_{\xi\xi'}^{\alpha}P_{\xi'}^{\alpha,i{\bf k}*} \]
where
\[ P_{\xi}^{\alpha,i{\bf k}} = \langle \psi_{i{\bf k}} | \beta_{\xi}^{\alpha} \rangle = \frac{1}{\Omega} \sum_{\bf G} \langle \psi_{i{\bf k}} | {\bf G+k} \rangle \langle {\bf G+k} | \beta_{\xi}^{\alpha} \rangle = \sum_{\bf G} \psi_i^{*}({\bf G+k}) \beta_{\xi}^{\alpha}({\bf G+k}) \]
where
\[ \beta_{\xi}^{\alpha}({\bf G+k}) = \frac{1}{\sqrt{\Omega}} \int e^{-i({\bf G+k}){\bf r}} \beta_{\xi}^{\alpha}({\bf r}) d{\bf r} = \frac{4\pi}{\sqrt{\Omega}} e^{-i{\bf (G+k) r_{\alpha}}}(-i)^{\ell} R_{\ell m}(\theta_{G+k}, \phi_{G+k}) \int \beta_{\ell}(r) j_{\ell}(|{\bf G+k}|r) r^2 dr \]
Contribution to stress tensor:
\[ \sigma_{\mu \nu}^{nl} = \frac{1}{\Omega} \frac{\partial E^{nl}}{\partial \varepsilon_{\mu \nu}} = \frac{1}{\Omega} \sum_{i{\bf k}} \sum_{\xi \xi'} \Bigg( \frac{\partial P_{\xi}^{\alpha,i{\bf k}}}{\partial \varepsilon_{\mu \nu}} D_{\xi\xi'}^{\alpha}P_{\xi'}^{\alpha,i{\bf k}*} + P_{\xi}^{\alpha,i{\bf k}} D_{\xi\xi'}^{\alpha} \frac{\partial P_{\xi'}^{\alpha,i{\bf k}*}}{\partial \varepsilon_{\mu \nu}} \Bigg) \]
We need to compute strain derivatives of \( P_{\xi}^{\alpha,i{\bf k}} \):
\[ \frac{\partial}{\partial \varepsilon_{\mu \nu}} P_{\xi}^{\alpha,i{\bf k}} = \sum_{\bf G} \psi_i^{*}({\bf G+k}) \frac{\partial}{\partial \varepsilon_{\mu \nu}} \beta_{\xi}^{\alpha}({\bf G+k}) \]
First way to take strain derivative of beta-projectors (here and below \( {\bf G+k} = {\bf q} \)):
\[ \frac{\partial}{\partial \varepsilon_{\mu \nu}} \beta_{\xi}^{\alpha}({\bf q}) = -\frac{4\pi}{2\sqrt{\Omega}} \delta_{\mu \nu} e^{-i{\bf q r_{\alpha}}}(-i)^{\ell} R_{\ell m}(\theta_q, \phi_q) \int \beta_{\ell}(r) j_{\ell}(q r) r^2 dr + \\ \frac{4\pi}{\sqrt{\Omega}} e^{-i{\bf q r_{\alpha}}}(-i)^{\ell} \frac{\partial R_{\ell m}(\theta_q, \phi_q)}{\partial \varepsilon_{\mu \nu}} \int \beta_{\ell}(r) j_{\ell}(q r) r^2 dr + \\ \frac{4\pi}{\sqrt{\Omega}} e^{-i{\bf q r_{\alpha}}}(-i)^{\ell} R_{\ell m}(\theta_q, \phi_q) \int \beta_{\ell}(r) \frac{\partial j_{\ell}(q r)}{\partial \varepsilon_{\mu \nu}} r^2 dr = \\ \frac{4\pi}{\sqrt{\Omega}} e^{-i{\bf q r_{\alpha}}}(-i)^{\ell} \Bigg[ \int \beta_{\ell}(r) j_{\ell}(q r) r^2 dr \Big(\frac{\partial R_{\ell m}(\theta_q, \phi_q)}{\partial \varepsilon_{\mu \nu}} - \frac{1}{2} R_{\ell m}(\theta_q, \phi_q) \delta_{\mu \nu}\Big) + R_{\ell m}(\theta_q, \phi_q) \int \beta_{\ell}(r) \frac{\partial j_{\ell}(q r)}{\partial \varepsilon_{\mu \nu}} r^2 dr \Bigg] \]
Strain derivative of the real spherical harmonics:
\[ \frac{\partial R_{\ell m}(\theta, \phi)}{\partial \varepsilon_{\mu \nu}} = \sum_{\tau} \frac{\partial R_{\ell m}(\theta, \phi)}{\partial q_{\tau}} \frac{\partial q_{\tau}}{\partial \varepsilon_{\mu \nu}} = -q_{\mu} \frac{\partial R_{\ell m}(\theta, \phi)}{\partial q_{\nu}} \]
For the derivatives of spherical harmonics over Cartesian components of vector please refer to the sht::dRlm_dr function.
Strain derivative of spherical Bessel function integral:
\[ \int \beta_{\ell}(r) \frac{\partial j_{\ell}(qr) }{\partial \varepsilon_{\mu \nu}} r^2 dr = \int \beta_{\ell}(r) \frac{\partial j_{\ell}(qr) }{\partial q} \frac{-q_{\mu} q_{\nu}}{q} r^2 dr \]
The second way to compute strain derivative of beta-projectors is trough Gaunt coefficients:
\[ \frac{\partial}{\partial \varepsilon_{\mu \nu}} \beta_{\xi}^{\alpha}({\bf q}) = -\frac{1}{2\sqrt{\Omega}} \delta_{\mu \nu} \int e^{-i{\bf q r}} \beta_{\xi}^{\alpha}({\bf r}) d{\bf r} + \frac{1}{\sqrt{\Omega}} \int i r_{\nu} q_{\mu} e^{-i{\bf q r}} \beta_{\xi}^{\alpha}({\bf r}) d{\bf r} \]
(second term comes from the strain derivative of \( e^{-i{\bf q r}} \)). Remembering that \( {\bf r} \) is proportional to p-like real spherical harmonics, we can rewrite the second part of beta-projector derivative as:
\[ \frac{1}{\sqrt{\Omega}} \int i r_{\nu} q_{\mu} e^{-i{\bf q r}} \beta_{\xi}^{\alpha}({\bf r}) d{\bf r} = \frac{1}{\sqrt{\Omega}} i q_{\mu} \int r \bar R_{1 \nu}(\theta, \phi) 4\pi \sum_{\ell_3 m_3} (-i)^{\ell_3} R_{\ell_3 m_3}(\theta_q, \phi_q) R_{\ell_3 m_3}(\theta, \phi) j_{\ell_3}(q r) \beta_{\ell_2}^{\alpha}(r) R_{\ell_2 m_2}(\theta, \phi) d{\bf r} = \\ \frac{4 \pi}{\sqrt{\Omega}} i q_{\mu} \sum_{\ell_3 m_3} (-i)^{\ell_3} R_{\ell_3 m_3}(\theta_q, \phi_q) \langle \bar R_{1\nu} | R_{\ell_3 m_3} | R_{\ell_2 m_2} \rangle \int j_{\ell_3}(q r) \beta_{\ell_2}^{\alpha}(r) r^3 dr \]
where
\[ \bar R_{1 x}(\theta, \phi) = -2 \sqrt{\frac{\pi}{3}} R_{11}(\theta, \phi) = \sin(\theta) \cos(\phi) \\ \bar R_{1 y}(\theta, \phi) = -2 \sqrt{\frac{\pi}{3}} R_{1-1}(\theta,\phi) = \sin(\theta) \sin(\phi) \\ \bar R_{1 z}(\theta, \phi) = 2 \sqrt{\frac{\pi}{3}} R_{10}(\theta, \phi) = \cos(\theta) \]
T | One of float, double, complex<float> or complex<double> types for generic or Gamma point case. |
Definition at line 37 of file stress.cpp.
r3::matrix< double > sirius::Stress::calc_stress_vloc | ( | ) |
Local potential contribution to stress.
Energy contribution from the local part of pseudopotential:
\[ E^{loc} = \int \rho({\bf r}) V^{loc}({\bf r}) d{\bf r} = \frac{1}{\Omega} \sum_{\bf G} \langle \rho | {\bf G} \rangle \langle {\bf G}| V^{loc} \rangle = \frac{1}{\Omega} \sum_{\bf G} \tilde \rho^{*}({\bf G}) \tilde V^{loc}({\bf G}) \]
where
\[ \tilde \rho({\bf G}) = \langle {\bf G} | \rho \rangle = \int e^{-i{\bf Gr}}\rho({\bf r}) d {\bf r} \]
and
\[ \tilde V^{loc}({\bf G}) = \langle {\bf G} | V^{loc} \rangle = \int e^{-i{\bf Gr}}V^{loc}({\bf r}) d {\bf r} \]
Using the expression for \( \tilde V^{loc}({\bf G}) \), the local contribution to the total energy is rewritten as
\[ E^{loc} = \frac{1}{\Omega} \sum_{\bf G} \tilde \rho^{*}({\bf G}) \sum_{\alpha} e^{-{\bf G\tau}_{\alpha}} 4 \pi \int V_{\alpha}^{loc}(r)\frac{\sin(Gr)}{Gr} r^2 dr = \frac{4\pi}{\Omega}\sum_{\bf G}\tilde \rho^{*}({\bf G})\sum_{\alpha} e^{-{\bf G\tau}_{\alpha}} \Bigg( \int \Big(V_{\alpha}(r) r + Z_{\alpha}^p {\rm erf}(r) \Big) \frac{\sin(Gr)}{G} dr - Z_{\alpha}^p \frac{e^{-\frac{G^2}{4}}}{G^2} \Bigg) \]
(see sirius::Potential::generate_local_potential for details).
Contribution to stress tensor:
\[ \sigma_{\mu \nu}^{loc} = \frac{1}{\Omega} \frac{\partial E^{loc}}{\partial \varepsilon_{\mu \nu}} = \frac{1}{\Omega} \frac{-1}{\Omega} \delta_{\mu \nu} \sum_{\bf G}\tilde \rho^{*}({\bf G}) \tilde V^{loc}({\bf G}) + \frac{4\pi}{\Omega^2} \sum_{\bf G}\tilde \rho^{*}({\bf G}) \sum_{\alpha} e^{-{\bf G\tau}_{\alpha}} \Bigg( \int \Big(V_{\alpha}(r) r + Z_{\alpha}^p {\rm erf}(r) \Big) \Big( \frac{r \cos (G r)}{G}-\frac{\sin (G r)}{G^2} \Big) \Big( -\frac{G_{\mu}G_{\nu}}{G} \Big) dr - Z_{\alpha}^p \Big(-\frac{e^{-\frac{G^2}{4}}}{2 G}-\frac{2 e^{-\frac{G^2}{4}}}{G^3} \Big) \Big( -\frac{G_{\mu}G_{\nu}}{G} \Big) \Bigg) = \\ -\delta_{\mu \nu} \sum_{\bf G}\rho^{*}({\bf G}) V^{loc}({\bf G}) + \sum_{\bf G} \rho^{*}({\bf G}) \Delta V^{loc}({\bf G}) G_{\mu}G_{\nu} \]
where \( \Delta V^{loc}({\bf G}) \) is built from the following radial integrals:
\[ \int \Big(V_{\alpha}(r) r + Z_{\alpha}^p {\rm erf}(r) \Big) \Big( \frac{\sin (G r)}{G^3} - \frac{r\cos (G r)}{G^2}\Big) dr - Z_{\alpha}^p \Big( \frac{e^{-\frac{G^2}{4}}}{2 G^2} + \frac{2 e^{-\frac{G^2}{4}}}{G^4}\Big) \]
Definition at line 706 of file stress.cpp.
|
inline |
Definition at line 279 of file stress.hpp.
r3::matrix< double > sirius::Stress::calc_stress_har | ( | ) |
Hartree energy contribution to stress.
Hartree energy:
\[ E^{H} = \frac{1}{2} \int_{\Omega} \rho({\bf r}) V^{H}({\bf r}) d{\bf r} = \frac{1}{2} \frac{1}{\Omega} \sum_{\bf G} \langle \rho | {\bf G} \rangle \langle {\bf G}| V^{H} \rangle = \frac{2 \pi}{\Omega} \sum_{\bf G} \frac{|\tilde \rho({\bf G})|^2}{G^2} \]
where
\[ \langle {\bf G} | \rho \rangle = \int e^{-i{\bf Gr}}\rho({\bf r}) d {\bf r} = \tilde \rho({\bf G}) \]
and
\[ \langle {\bf G} | V^{H} \rangle = \int e^{-i{\bf Gr}}V^{H}({\bf r}) d {\bf r} = \frac{4\pi}{G^2} \tilde \rho({\bf G}) \]
Hartree energy contribution to stress tensor:
\[ \sigma_{\mu \nu}^{H} = \frac{1}{\Omega} \frac{\partial E^{H}}{\partial \varepsilon_{\mu \nu}} = \frac{1}{\Omega} 2\pi \Big( \big( \frac{\partial}{\partial \varepsilon_{\mu \nu}} \frac{1}{\Omega} \big) \sum_{{\bf G}} \frac{|\tilde \rho({\bf G})|^2}{G^2} + \frac{1}{\Omega} \sum_{{\bf G}} |\tilde \rho({\bf G})|^2 \frac{\partial}{\partial \varepsilon_{\mu \nu}} \frac{1}{G^2} \Big) = \\ \frac{1}{\Omega} 2\pi \Big( -\frac{1}{\Omega} \delta_{\mu \nu} \sum_{{\bf G}} \frac{|\tilde \rho({\bf G})|^2}{G^2} + \frac{1}{\Omega} \sum_{\bf G} |\tilde \rho({\bf G})|^2 \sum_{\tau} \frac{-2 G_{\tau}}{G^4} \frac{\partial G_{\tau}}{\partial \varepsilon_{\mu \nu}} \Big) = \\ 2\pi \sum_{\bf G} \frac{|\rho({\bf G})|^2}{G^2} \Big( -\delta_{\mu \nu} + \frac{2}{G^2} G_{\nu} G_{\mu} \Big) \]
Definition at line 613 of file stress.cpp.
|
inline |
Definition at line 315 of file stress.hpp.
r3::matrix< double > sirius::Stress::calc_stress_ewald | ( | ) |
Ewald energy contribution to stress.
Ewald energy:
\[ E^{ion-ion} = \frac{1}{2} \sideset{}{'} \sum_{\alpha \beta {\bf T}} Z_{\alpha} Z_{\beta} \frac{{\rm erfc}(\sqrt{\lambda} |{\bf r}_{\alpha} - {\bf r}_{\beta} + {\bf T}|)}{|{\bf r}_{\alpha} - {\bf r}_{\beta} + {\bf T}|} + \frac{2 \pi}{\Omega} \sum_{{\bf G}} \frac{e^{-\frac{G^2}{4 \lambda}}}{G^2} \Big| \sum_{\alpha} Z_{\alpha} e^{-i{\bf r}_{\alpha}{\bf G}} \Big|^2 - \sum_{\alpha} Z_{\alpha}^2 \sqrt{\frac{\lambda}{\pi}} - \frac{2\pi}{\Omega}\frac{N_{el}^2}{4 \lambda} \]
(check sirius::DFT_ground_state::ewald_energy for details).
Contribution to stress tensor:
\[ \sigma_{\mu \nu}^{ion-ion} = \frac{1}{\Omega} \frac{\partial E^{ion-ion}}{\partial \varepsilon_{\mu \nu}} \]
Derivative of the first part:
\[ \frac{1}{\Omega}\frac{\partial}{\partial \varepsilon_{\mu \nu}} \frac{1}{2} \sideset{}{'} \sum_{\alpha \beta {\bf T}} Z_{\alpha} Z_{\beta} \frac{{\rm erfc}(\sqrt{\lambda} |{\bf r}_{\alpha} - {\bf r}_{\beta} + {\bf T}|)}{|{\bf r}_{\alpha} - {\bf r}_{\beta} + {\bf T}|} = \frac{1}{2\Omega} \sideset{}{'} \sum_{\alpha \beta {\bf T}} Z_{\alpha} Z_{\beta} \Big( -2e^{-\lambda |{\bf r'}|^2} \sqrt{\frac{\lambda}{\pi}} \frac{1}{|{\bf r'}|^2} - {\rm erfc}(\sqrt{\lambda} |{\bf r'}|) \frac{1}{|{\bf r'}|^3} \Big) r'_{\mu} r'_{\nu} \]
where \( {\bf r'} = {\bf r}_{\alpha} - {\bf r}_{\beta} + {\bf T} \).
Derivative of the second part:
\[ \frac{1}{\Omega}\frac{\partial}{\partial \varepsilon_{\mu \nu}} \frac{2\pi}{\Omega} \sum_{{\bf G}} \frac{e^{-\frac{G^2}{4 \lambda}}}{G^2} \Big| \sum_{\alpha} Z_{\alpha} e^{-i{\bf r}_{\alpha}{\bf G}} \Big|^2 = -\frac{2\pi}{\Omega^2} \delta_{\mu \nu} \sum_{{\bf G}} \frac{e^{-\frac{G^2}{4 \lambda}}}{G^2} \Big| \sum_{\alpha} Z_{\alpha} e^{-i{\bf r}_{\alpha}{\bf G}} \Big|^2 + \frac{2\pi}{\Omega^2} \sum_{\bf G} G_{\mu} G_{\nu} \frac{e^{-\frac{G^2}{4 \lambda}}}{G^2} 2 \frac{\frac{G^2}{4\lambda} + 1}{G^2} \Big| \sum_{\alpha} Z_{\alpha} e^{-i{\bf r}_{\alpha}{\bf G}} \Big|^2 \]
Derivative of the fourth part:
\[ -\frac{1}{\Omega}\frac{\partial}{\partial \varepsilon_{\mu \nu}} \frac{2\pi}{\Omega}\frac{N_{el}^2}{4 \lambda} = \frac{2\pi}{\Omega^2}\frac{N_{el}^2}{4 \lambda} \delta_{\mu \nu} \]
Definition at line 490 of file stress.cpp.
|
inline |
Definition at line 362 of file stress.hpp.
void sirius::Stress::calc_stress_kin_aux |
Kinetic energy contribution to stress.
Kinetic energy:
\[ E^{kin} = \sum_{{\bf k}} w_{\bf k} \sum_j f_j \frac{1}{2} |{\bf G+k}|^2 |\psi_j({\bf G + k})|^2 \]
Contribution to the stress tensor
\[ \sigma_{\mu \nu}^{kin} = \frac{1}{\Omega} \frac{\partial E^{kin}}{\partial \varepsilon_{\mu \nu}} = \frac{1}{\Omega} \sum_{{\bf k}} w_{\bf k} \sum_j f_j \frac{1}{2} 2 |{\bf G+k}| \Big( -\frac{1}{|{\bf G+k}|} (G+k)_{\mu} (G+k)_{\nu} \Big) |\psi_j({\bf G + k})|^2 =\\ -\frac{1}{\Omega} \sum_{{\bf k}} w_{\bf k} (G+k)_{\mu} (G+k)_{\nu} \sum_j f_j |\psi_j({\bf G + k})|^2 \]
Definition at line 649 of file stress.cpp.
r3::matrix< double > sirius::Stress::calc_stress_kin | ( | ) |
Definition at line 692 of file stress.cpp.
|
inline |
Definition at line 385 of file stress.hpp.
r3::matrix< double > sirius::Stress::calc_stress_nonloc | ( | ) |
Definition at line 758 of file stress.cpp.
|
inline |
Definition at line 392 of file stress.hpp.
r3::matrix< double > sirius::Stress::calc_stress_us | ( | ) |
Contribution to the stress tensor from the augmentation operator.
Total energy in ultrasoft pseudopotential contains this term:
\[ \int V^{eff}({\bf r})\rho^{aug}({\bf r})d{\bf r} = \sum_{\alpha} \sum_{\xi \xi'} n_{\xi \xi'}^{\alpha} \int V^{eff}({\bf r}) Q_{\xi \xi'}^{\alpha}({\bf r}) d{\bf r} \]
The derivatives of beta-projectors (hidden in the desnity matrix expression) are taken into account in sirius::Stress::calc_stress_nonloc. Here we need to compute the remaining contribution from the \( Q_{\xi \xi'}({\bf r}) \) itself. We are interested in the integral:
\[ \int V^{eff}({\bf r}) Q_{\xi \xi'}^{\alpha}({\bf r}) d{\bf r} = \sum_{\bf G} V^{eff}({\bf G}) \tilde Q_{\xi \xi'}^{\alpha}({\bf G}) \]
where
\[ \tilde Q_{\xi \xi'}^{\alpha}({\bf G}) = \int e^{-i{\bf Gr}} Q_{\xi \xi'}^{\alpha}({\bf r}) d{\bf r} = 4\pi \sum_{\ell m} (-i)^{\ell} R_{\ell m}(\hat{\bf G}) \langle R_{\ell_{\xi} m_{\xi}} | R_{\ell m} | R_{\ell_{\xi'} m_{\xi'}} \rangle \int Q_{\ell_{\xi} \ell_{\xi'}}^{\ell}(r) j_{\ell}(Gr) r^2 dr \]
Strain derivative of \( \tilde Q_{\xi \xi'}^{\alpha}({\bf G}) \) is:
\[ \frac{\partial}{\partial \varepsilon_{\mu \nu}} \tilde Q_{\xi \xi'}^{\alpha}({\bf G}) = 4\pi \sum_{\ell m} (-i)^{\ell} \langle R_{\ell_{\xi} m_{\xi}} | R_{\ell m} | R_{\ell_{\xi'} m_{\xi'}} \rangle \Big( \frac{\partial R_{\ell m}(\hat{\bf G})}{\partial \varepsilon_{\mu \nu}} \int Q_{\ell_{\xi} \ell_{\xi'}}^{\ell}(r) j_{\ell}(Gr) r^2 dr + R_{\ell m}(\hat{\bf G}) \int Q_{\ell_{\xi} \ell_{\xi'}}^{\ell}(r) \frac{\partial j_{\ell}(Gr)}{\partial \varepsilon_{\mu \nu}} r^2 dr \Big) \]
For strain derivatives of spherical harmonics and Bessel functions see sirius::Stress::calc_stress_nonloc. We can pull the common multiplier \( -G_{\mu} / G \) from both terms and arrive to the following expression:
\[ \frac{\partial}{\partial \varepsilon_{\mu \nu}} \tilde Q_{\xi \xi'}^{\alpha}({\bf G}) = -\frac{G_{\mu}}{G} 4\pi \sum_{\ell m} (-i)^{\ell} \langle R_{\ell_{\xi} m_{\xi}} | R_{\ell m} | R_{\ell_{\xi'} m_{\xi'}} \rangle \Big( \big(\nabla_{G} R_{\ell m}(\hat{\bf G})\big)_{\nu} \int Q_{\ell_{\xi} \ell_{\xi'}}^{\ell}(r) j_{\ell}(Gr) r^2 dr + R_{\ell m}(\hat{\bf G}) \int Q_{\ell_{\xi} \ell_{\xi'}}^{\ell}(r) \frac{\partial j_{\ell}(Gr)}{\partial G} G_{\nu} r^2 dr \Big) \]
Definition at line 371 of file stress.cpp.
|
inline |
Definition at line 441 of file stress.hpp.
|
inline |
Definition at line 446 of file stress.hpp.
r3::matrix< double > sirius::Stress::calc_stress_xc | ( | ) |
XC contribution to stress.
XC contribution has the following expression:
\[ \frac{\partial E_{xc}}{\partial \varepsilon_{\mu \nu}} = \delta_{\mu \nu} \int \Big( \epsilon_{xc}({\bf r}) - v_{xc}({\bf r}) \Big) \rho({\bf r})d{\bf r} - \int \frac{\partial \epsilon_{xc} \big( \rho({\bf r}), \nabla \rho({\bf r})\big) }{\nabla_{\mu} \rho({\bf r})} \nabla_{\nu}\rho({\bf r}) d{\bf r} \]
Definition at line 284 of file stress.cpp.
|
inline |
Definition at line 461 of file stress.hpp.
r3::matrix< double > sirius::Stress::calc_stress_core | ( | ) |
Non-linear core correction to stress tensor.
Definition at line 217 of file stress.cpp.
|
inline |
Definition at line 469 of file stress.hpp.
r3::matrix< double > sirius::Stress::calc_stress_hubbard | ( | ) |
Definition at line 114 of file stress.cpp.
|
inline |
Definition at line 476 of file stress.hpp.
r3::matrix< double > sirius::Stress::calc_stress_total | ( | ) |
Definition at line 87 of file stress.cpp.
|
inline |
Definition at line 483 of file stress.hpp.
void sirius::Stress::print_info | ( | std::ostream & | out__, |
int | verbosity__ | ||
) | const |
Definition at line 561 of file stress.cpp.
|
private |
Definition at line 99 of file stress.hpp.
|
private |
Definition at line 101 of file stress.hpp.
|
private |
Definition at line 103 of file stress.hpp.
|
private |
Definition at line 105 of file stress.hpp.
|
private |
Definition at line 107 of file stress.hpp.
|
private |
Definition at line 109 of file stress.hpp.
|
private |
Definition at line 111 of file stress.hpp.
|
private |
Definition at line 113 of file stress.hpp.
|
private |
Definition at line 115 of file stress.hpp.
|
private |
Definition at line 117 of file stress.hpp.
|
private |
Definition at line 119 of file stress.hpp.
|
private |
Definition at line 121 of file stress.hpp.
|
private |
Definition at line 123 of file stress.hpp.
|
private |
Definition at line 125 of file stress.hpp.