Octopus
|
This module implements the calculation of the stress tensor. More...
This module implements the calculation of the stress tensor.
Functions/Subroutines | |
subroutine, public | stress_calculate (namespace, gr, hm, st, ions, ks, ext_partners) |
This computes the total stress on the lattice. More... | |
subroutine | stress_from_kinetic (gr, space, hm, st, symm, rcell_volume, stress_kin) |
Computes the contribution to the stress tensor from the kinetic energy. More... | |
subroutine | stress_from_hartree (gr, space, volume, vh, grad_vh, ehartree, stress_Hartree) |
Computes the contribution to the stress tensor from the Hartree energy. More... | |
subroutine | stress_from_xc (energy, rcell_volume, periodic_dim, stress_xc) |
Computes the contribution to the stress tensor from the xc energy. More... | |
subroutine | stress_from_xc_nlcc (rcell_volume, gr, st, ions, vxc, stress_xc_nlcc) |
Computes the NLCC contribution to the stress tensor from the xc energy. More... | |
subroutine | stress_from_pseudo_nonloc (gr, st, hm, ions, stress_ps_nl) |
Computes the contribution to the stress tensor from the nonlocal part of the pseudopotentials. More... | |
subroutine | stress_from_pseudo_local (gr, st, hm, ions, rho_total, vh, grad_vh, stress_ps_local) |
Computes the contribution from the local part of the pseudopotential. More... | |
subroutine | epot_local_pseudopotential_sr (mesh, ions, iatom, vpsl, rvpsl) |
subroutine | stress_from_hubbard (namespace, gr, st, hm, space, rcell_volume, stress_hubbard) |
Computes the contribution to the stress tensor from the Hubbard energy. More... | |
subroutine, public | output_stress (iunit, space_dim, stress_tensors, all_terms) |
subroutine, public | output_pressure (iunit, space_dim, total_stress_tensor) |
subroutine | print_stress_tensor (ounit, space_dim, tensor) |
subroutine, public stress_oct_m::stress_calculate | ( | type(namespace_t), intent(in) | namespace, |
type(grid_t), intent(inout) | gr, | ||
type(hamiltonian_elec_t), intent(inout) | hm, | ||
type(states_elec_t), intent(inout) | st, | ||
type(ions_t), intent(inout) | ions, | ||
type(v_ks_t), intent(in) | ks, | ||
type(partner_list_t), intent(in) | ext_partners | ||
) |
This computes the total stress on the lattice.
[in,out] | gr | grid |
[in,out] | hm | the Hamiltonian |
[in,out] | st | the electronic states |
[in,out] | ions | geometry |
[in] | ks | the Kohn-Sham system |
[in] | ext_partners | external interaction partners |
Definition at line 183 of file stress.F90.
|
private |
Computes the contribution to the stress tensor from the kinetic energy.
We use the real space formula from Sharma and Suryanarayana On the calculation of the stress tensor in real-space Kohn-Sham density functional theory J. Chem. Phys. 149, 194104 (2018)
More precisely, this routines computes
\[ \sigma_{ij} = \frac{1}{V}\sum_n\sum_k^{BZ} w_kf_{n,k}\int d^3r \partial_i \psi^*_{n,k}(r) \partial_j \psi_{n,k}(r)\, \]
where \(V\) is the cell volume, \( w_k\) is the weight of the k-point k, \( f_{n,k}\) is the occupation number of the band with a k-point index k, and \( \psi_{n,k}\) is the corresponding Bloch state.
Definition at line 389 of file stress.F90.
|
private |
Computes the contribution to the stress tensor from the Hartree energy.
We use the real space formula from Sharma and Suryanarayana On the calculation of the stress tensor in real-space Kohn-Sham density functional theory J. Chem. Phys. 149, 194104 (2018)
More precisely, this routines computes
\[ \sigma_{ij} = \frac{1}{4\pi V}\int d^3r \partial_i v_{\rm H}(r) \partial_j v_{\rm H}(r) - \delta_{ij} \frac{1}{V}E_{\rm H}\, \]
where \(V\) is the cell volume, \( v_{\rm H}(r) \) is the Hartree potential, and \( E_{\rm H} \) is the Hartree energy defined as
\[ E_{\rm H} = \frac{1}{2} \int d^3r n(r) v_{\rm H}(r) \,. \]
[in] | vh | Hartree potential |
[in] | grad_vh | Gradient of the Hartree potential |
[in] | ehartree | Hartree U = (1/2)*Int [n v_Hartree] |
Definition at line 478 of file stress.F90.
|
private |
Computes the contribution to the stress tensor from the xc energy.
This routines computes the xc stress tensor assuming an LDA functional
\[ \sigma_{ij} = \frac{1}{V}\delta_{ij}(-\int d^3r n(r) v_{\rm xc}(r) + E_{\rm xc})\, \]
where \(V\) is the cell volume, \(n(r) \) is the electronic density, \( v_{\rm xc}(r) \) is the exchange-correlation potential, and \( E_{\rm xc} \) is the exchange-correlation energy.
The GGA extra term is computed alongside with the calculation of the xc potential, in the routine xc_get_vxc.
Definition at line 527 of file stress.F90.
|
private |
Computes the NLCC contribution to the stress tensor from the xc energy.
The nonlinear core correction term is given by
\[ \sigma_{ij}^{\rm NLCC} = \frac{1}{V}\int d^3r v_{xc}(r) \frac{\partial \rho_{\rm NLCC}(\epsilon r)}{\partial \epsilon_{ij}}\Bigg|_{\epsilon=I}\,. \]
Definition at line 556 of file stress.F90.
|
private |
Computes the contribution to the stress tensor from the nonlocal part of the pseudopotentials.
More precisely, this routines computes
\[ \sigma_{ij} = \frac{2}{V}\sum_n\sum_k^{BZ} \sum_I w_kf_{n,k} \langle \partial_i \psi_{n,k}| (r_j-R_j) \hat{V}_{\rm NL,I}| \psi_{n,k}\rangle + \delta_{ij} \frac{1}{V}E_{\rm NL} \]
where \(V\) is the cell volume, \( w_k\) is the weight of the k-point k, \( f_{n,k}\) is the occupation number of the band with a k-point index k, \( \psi_{n,k}\) is the corresponding Bloch state, and \( \hat{V}_{\rm NL, I}\) is the pseudopotential non-local operator from atom I, centered on the position \(R_I\). \( E_{\rm NL} \) is the nonlocal energy from the nonlocal pseudopotential.
See Sharma and Suryanarayana, On the calculation of the stress tensor in real-space Kohn-Sham density functional theory, J. Chem. Phys. 149, 194104 (2018) for more details
[in] | gr | grid |
Definition at line 629 of file stress.F90.
|
private |
Computes the contribution from the local part of the pseudopotential.
We use a real space formulation, which computes two parts. One short range (SR)
\[ \sigma_{ij} =\frac{1}{V} \sum_I\int d^3r [\partial_i \rho(r)] (x_j-R_{I,j}) v_{\rm SR, I}(r) + \delta_{ij} \frac{1}{V}E_{\rm loc, SR}\, \]
and a long range part
\[ \sigma_{ij} = \frac{1}{V}\int d\vec{x} \sum_I n_{\rm LR}^{I}(|\vec{x}-\vec{R}_I|) (x_j-R_{I,j}) [\partial_i v_H(\vec{x})] + \delta_{ij} \frac{1}{V} E_{\rm loc, LR} -\frac{2}{4\pi V}\int d\vec{r} [\partial_i v_{\rm H}(\vec{r})][\partial_j v_{\rm LR}(\vec{r})] \]
where \(V\) is the cell volume, \( v_{\rm H}(r) \) is the Hartree potential, \(\vec{R}_I\) refers to the position of the atom \(I\), and \(n_{\rm LR}^{I}\) is the long-range density associated with the long-range potential \(v_{\rm LR}^I\).
[in] | gr | grid |
[in] | vh | Hartree potential |
[in] | grad_vh | Gradient of the Hartree potential |
Definition at line 735 of file stress.F90.
|
private |
Definition at line 874 of file stress.F90.
|
private |
Computes the contribution to the stress tensor from the Hubbard energy.
More precisely, this routine computes
\[ \sigma^{U}_{\alpha\beta} = \delta_{\alpha\beta} E_{\rm Hubbard} - \frac{2}{V} \sum_{n} \sum_{\mathbf{k}}^{BZ} w_k \Re\Bigg\{ \langle [\partial_\alpha \psi^*_{n\mathbf{k}}(\mathbf{r})] | (\mathbf{r}-\mathbf{R}_J)_\beta \hat{V}_U n_\mathbf{k}\rangle\Bigg\} \]
where \(V\) is the cell volume, \( w_k\) is the weight of the k-point k, \( \psi_{n,\mathbf{k}}\) is the corresponding Bloch state, and \( \hat{V}_{U}\) is the Hubbard potential and \( E_{\rm Hubbard} \) is the Hubbard energy.
[in] | gr | grid |
Definition at line 938 of file stress.F90.
subroutine, public stress_oct_m::output_stress | ( | integer, intent(in) | iunit, |
integer, intent(in) | space_dim, | ||
type(stress_t), intent(in) | stress_tensors, | ||
logical, intent(in), optional | all_terms | ||
) |
[in] | all_terms | if yes, writes each contributing term separately |
Definition at line 1032 of file stress.F90.
subroutine, public stress_oct_m::output_pressure | ( | integer, intent(in) | iunit, |
integer, intent(in) | space_dim, | ||
real(real64), dimension(3,3), intent(in) | total_stress_tensor | ||
) |
Definition at line 1098 of file stress.F90.
|
private |
Definition at line 1127 of file stress.F90.