Octopus
|
Data Types | |
type | ion_interaction_t |
Functions/Subroutines | |
subroutine, public | ion_interaction_init (this, namespace, space, natoms) |
subroutine, public | ion_interaction_init_parallelization (this, natoms, mc) |
subroutine, public | ion_interaction_end (this) |
subroutine, public | ion_interaction_calculate (this, space, latt, atom, natoms, pos, lsize, energy, force, energy_components, force_components) |
Top level routine for computing electrostatic energies and forces between ions. More... | |
real(real64) function | jellium_slab_energy_periodic (space, atom, lsize) |
Electrostatic energy of a periodic jellium slab. More... | |
real(real64) function | jellium_self_energy_finite (dist, latt, atom, lsize) |
Electrostatic self-interaction for jellium instances, with orthogonal cells. More... | |
subroutine | ion_interaction_finite (dist, space, atom, pos, lsize, energy, force) |
Electrostatic Ewald energy and forces for finite systems. More... | |
subroutine | ion_interaction_periodic (this, space, latt, atom, natoms, pos, energy, force, energy_components, force_components) |
Total Ewald electrostatic energy and forces, for 1D, 2D and 3D systems. More... | |
subroutine | ewald_short (dist, space, latt, atom, pos, alpha, ereal, force) |
Short range component of the Ewald electrostatic energy and force. More... | |
subroutine | ewald_self_interaction (dist, atom, alpha, eself, charge) |
@ brief Ewald self-interaction energy More... | |
subroutine | ewald_long_3d (this, space, latt, atom, natoms, pos, efourier, force, charge) |
Computes the long-range part of the 3D Ewald summation. More... | |
subroutine | ewald_long_2d (this, space, latt, atom, natoms, pos, efourier, force) |
In-Chul Yeh and Max L. Berkowitz, J. Chem. Phys. 111, 3155 (1999). More... | |
subroutine | pseudopotential_correction_3d (dist, latt, atom, charge, epseudo) |
G=0 component of Ewald energy arising from the pseudopotentials, for 3D systems. More... | |
subroutine, public | ion_interaction_stress (this, space, latt, atom, natoms, pos, stress_ii) |
Computes the contribution to the stress tensor the ion-ion energy. More... | |
subroutine | ion_interaction_stress_short (this, space, latt, atom, natoms, pos, stress_short) |
Computes the short-range contribution to the stress tensor the ion-ion energy. More... | |
subroutine | ewald_3d_stress (this, space, latt, atom, natoms, pos, stress_Ewald) |
Computes the contribution to the stress tensor from the 3D Ewald sum. More... | |
subroutine | ewald_2d_stress (this, space, latt, atom, natoms, pos, stress_Ewald) |
Computes the contribution to the stress tensor from the 2D Ewald sum. More... | |
real(real64) function | screening_function_2d (alpha, z_ij, gg_abs, erfc) |
Auxiliary function for the Ewald 2D stress. More... | |
subroutine, public | ion_interaction_test (space, latt, atom, natoms, pos, lsize, namespace, mc) |
Variables | |
integer, parameter | ion_component_real = 1 |
integer, parameter | ion_component_self = 2 |
integer, parameter | ion_component_fourier = 3 |
integer, parameter | ion_num_components = 3 |
subroutine, public ion_interaction_oct_m::ion_interaction_init | ( | type(ion_interaction_t), intent(out) | this, |
type(namespace_t), intent(in) | namespace, | ||
class(space_t), intent(in) | space, | ||
integer, intent(in) | natoms | ||
) |
Definition at line 162 of file ion_interaction.F90.
subroutine, public ion_interaction_oct_m::ion_interaction_init_parallelization | ( | type(ion_interaction_t), intent(inout) | this, |
integer, intent(in) | natoms, | ||
type(multicomm_t), intent(in) | mc | ||
) |
Definition at line 194 of file ion_interaction.F90.
subroutine, public ion_interaction_oct_m::ion_interaction_end | ( | type(ion_interaction_t), intent(inout) | this | ) |
Definition at line 212 of file ion_interaction.F90.
subroutine, public ion_interaction_oct_m::ion_interaction_calculate | ( | type(ion_interaction_t), intent(inout) | this, |
class(space_t), intent(in) | space, | ||
type(lattice_vectors_t), intent(in) | latt, | ||
type(atom_t), dimension(:), intent(in) | atom, | ||
integer, intent(in) | natoms, | ||
real(real64), dimension(1:space%dim,1:natoms), intent(in) | pos, | ||
real(real64), dimension(:), intent(in) | lsize, | ||
real(real64), intent(out) | energy, | ||
real(real64), dimension(:, :), intent(out) | force, | ||
real(real64), dimension(:), intent(out), optional | energy_components, | ||
real(real64), dimension(:, :, :), intent(out), optional | force_components | ||
) |
Top level routine for computing electrostatic energies and forces between ions.
For additional details about this routine, see http:
[in,out] | this | Ion interaction instance |
[in] | space | System dimensions and boundary conditions |
[in] | latt | Crystal lattice |
[in] | atom | Atoms |
[in] | natoms | Number of atoms == size(atoms) |
[in] | pos | Atomic positions |
[in] | lsize | Box half-lengths |
[out] | energy | Total ion-ion electrostatic energy |
[out] | force | Total force on each ion |
[out] | energy_components | Energy contributions |
[out] | force_components | Force contributions |
Definition at line 228 of file ion_interaction.F90.
|
private |
Electrostatic energy of a periodic jellium slab.
From the inner potential been \( v(z) = 4\pi \rho_0 ( z^2/2 + L^2/8)\) we get the energy as \(U = 1/2 \int \rho(r) v(r) d^3r = \pi / 3 \rho^2 A L^3\)
[in] | space | System dimensions and boundary conditions |
[in] | atom | Atoms |
[in] | lsize | Box half-lengths |
Definition at line 282 of file ion_interaction.F90.
|
private |
Electrostatic self-interaction for jellium instances, with orthogonal cells.
If no atoms are of species SPECIES_JELLIUM or SPECIES_JELLIUM_SLAB, the returned energy is zero.
[in] | dist | Atom MPI distribution |
[in] | latt | Crystal lattice |
[in] | atom | Atoms |
[in] | lsize | Box half-lengths |
Definition at line 318 of file ion_interaction.F90.
|
private |
Electrostatic Ewald energy and forces for finite systems.
[in] | dist | Atom distribution |
[in] | space | System dimensions |
[in] | atom | Atoms |
[in] | pos | Ion positions |
[in] | lsize | Box half-lengths |
[out] | energy | Total electrostatic energy |
[out] | force | Forces |
Definition at line 361 of file ion_interaction.F90.
|
private |
Total Ewald electrostatic energy and forces, for 1D, 2D and 3D systems.
[in] | this | Ion interaction instance |
[in] | space | System dimensions and boundary conditions |
[in] | latt | Crystal lattice |
[in] | atom | Atoms |
[in] | natoms | Number of atoms in system |
[in] | pos | Atomic positions |
[out] | energy | Total ion-ion electrostatic energy |
[out] | force | Total force on each ion |
[out] | energy_components | Ion-ion energy contributions |
[out] | force_components | Force contributions |
Definition at line 412 of file ion_interaction.F90.
|
private |
Short range component of the Ewald electrostatic energy and force.
Computes the energy:
\[ E_{\text{short}} = \frac{1}{2} \sum_{i \neq j}^{N_{atom}} \sum_{\mathbf{T}} \frac{Z_i Z_j \operatorname{erfc}\left(\sqrt{\alpha} |\mathbf{r}_{i j} + \mathbf{T}|\right)}{|\mathbf{r}_{i j} + \mathbf{T}|}, \]
and the force:
\[ \mathbf{f}_j = \sum_i \sum_{\mathbf{T}} -Z_i Z_j (\mathbf{r}_{ij} + \mathbf{T}) \frac{\operatorname{erfc}}{|\mathbf{r}_{ij} + \mathbf{T}|} + \frac{2\alpha}{\sqrt{\pi}}\frac{e^{-\alpha |\mathbf{r}_{ij} + \mathbf{T}|^2}}{|\mathbf{r}_{ij} + \mathbf{T}|^2} \]
The middle loop over atoms only iterates over the upper triangle, resulting in a large load imbalance for MPI, however it should still be an improvement w.r.t. iterating over all atoms. OMP cannot be used to collapse non-rectangular loops.
One could alternatively collapse the loops over the whole matrix and distribute a single index, so each process would get its own set of (ia, ja, icopy). The individual indices could be recovered algebraically.
[in] | dist | Atomic distribution |
[in] | space | System dims and boundary conditions |
[in] | latt | Lattice |
[in] | atom | Atom properties |
[in] | pos | Ion positions |
[in] | alpha | Broadening parameter, controlling range separation |
[out] | ereal | Realspace contribution to Ewald energy |
[in,out] | force | Forces. |
Definition at line 500 of file ion_interaction.F90.
|
private |
@ brief Ewald self-interaction energy
The force exerted on an atom by itself is by definition zero. Also returns the total ion charge.
[in] | dist | Atomic distribution |
[in] | atom | Atom properties |
[in] | alpha | Broadening parameter, controlling range separation |
[out] | eself | Self-interaction energy |
[out] | charge | Total ionic charge |
Definition at line 582 of file ion_interaction.F90.
|
private |
Computes the long-range part of the 3D Ewald summation.
[in] | this | Ion interaction instance |
[in] | space | System dimensions and boundary conditions |
[in] | latt | Crystal lattice |
[in] | atom | Atoms |
[in] | natoms | Number of atoms in system |
[in] | pos | (spacedim, natoms) |
[in,out] | efourier | Long-range part of total ion-ion energy |
[in,out] | force | (spacedim, natoms) |
[in] | charge | Total ionic charge |
Definition at line 610 of file ion_interaction.F90.
|
private |
In-Chul Yeh and Max L. Berkowitz, J. Chem. Phys. 111, 3155 (1999).
[in,out] | force | (spacedim, natoms) |
Definition at line 696 of file ion_interaction.F90.
|
private |
G=0 component of Ewald energy arising from the pseudopotentials, for 3D systems.
See J. Ihm, A. Zunger, M.L. Cohen, J. Phys. C 12, 4409 (1979) This is the term \alpha_1 or Eq. 12.22 in the book of R. Martin.
[in] | dist | Atom MPI distribution |
[in] | latt | Lattice |
[in] | atom | Atoms |
[out] | epseudo | Ewald energy correction |
charge | Total ionic charge |
Definition at line 901 of file ion_interaction.F90.
subroutine, public ion_interaction_oct_m::ion_interaction_stress | ( | type(ion_interaction_t), intent(inout) | this, |
class(space_t), intent(in) | space, | ||
type(lattice_vectors_t), intent(in) | latt, | ||
type(atom_t), dimension(:), intent(in) | atom, | ||
integer, intent(in) | natoms, | ||
real(real64), dimension(:,:), intent(in) | pos, | ||
real(real64), dimension(space%dim, space%dim), intent(out) | stress_ii | ||
) |
Computes the contribution to the stress tensor the ion-ion energy.
Definition at line 929 of file ion_interaction.F90.
|
private |
Computes the short-range contribution to the stress tensor the ion-ion energy.
The formula that is implemented here corresponds to the short-range part of the expression B1 of Nielsen and Martin, Stresses in semiconductors: Ab initio calculations on Si, Ge, and GaAs PRB 32, 3792 (1985).
\[ \sigma_{\alpha\beta}^{\rm Ewald-SR} = \frac{\sqrt{\epsilon}}{2\Omega} \sum_{\mathbf{\tau\sigma T}} Z_\tau Z_{\sigma} H(\sqrt{\epsilon}D) \frac{D_\alpha D_{\beta}}{D^2}\Bigg|_{\mathbf{D=x_{\sigma}-x_{\tau}+T\neq0}}\,, \]
where the function \(H(x)\) is
\[ H(x) = \frac{\partial[{\rm erfc}(x)]}{\partial x} - \frac{{\rm erfc}(x)}{x} \]
Definition at line 982 of file ion_interaction.F90.
|
private |
Computes the contribution to the stress tensor from the 3D Ewald sum.
The formula that is implemented here correspond to the long-range part of the expression B1 of Nielsen and Martin, Stresses in semiconductors: Ab initio calculations on Si, Ge, and GaAs PRB 32, 3792 (1985).
\[ \sigma_{\alpha\beta}^{\rm Ewald-LR} = \frac{\pi}{2\Omega^2\epsilon} \sum_{\mathbf{G\neq0}} \frac{e^{-G^2/4\epsilon}}{G^2/4\epsilon}\left|\sum_\tau Z_\tau e^{i\mathbf{G}\cdot\mathbf{x}_\tau}\right|^2 \Big(2\frac{G_\alpha G_\beta}{G^2}(G^2/4\epsilon + 1) - \delta_{\alpha\beta}\Big) \nonumber\ + \frac{\pi}{2\Omega^2\epsilon}\left(\sum_{\tau} Z_\tau\right)^2\delta_{\alpha\beta}\,, \]
Definition at line 1067 of file ion_interaction.F90.
|
private |
Computes the contribution to the stress tensor from the 2D Ewald sum.
The formula that is implemented here is derived from the corresponding energy
\[ \frac{\partial \gamma_{\rm Ewald, LR-2D}}{\partial \epsilon_{\alpha\beta}}\Bigg|_{\vec{\epsilon}=\vec{I}} = - \delta_{\alpha\beta} \gamma_{\rm Ewald, LR-2D} + \frac{\pi}{4A\alpha}\sum_{ij}Z_i Z_{j}\sum_{G\neq 0}\cos(\vec{G}\cdot\vec{r}_{ij})\frac{G_\alpha G_\beta}{G^2}\nonumber\ \Bigg( e^{G z_{ij}}\Bigg[ 2\alpha\Big(\frac{1}{G}-z_{ij}){\rm erfc}(\alpha z_{ij} + \frac{G}{2\alpha}) -{\rm erfc`}(\alpha z_{ij} + \frac{G}{2\alpha})\Bigg]\nonumber\ + e^{-G z_{ij}}\Bigg[ 2\alpha\Big(\frac{1}{G}+z_{ij}){\rm erfc}(-\alpha z_{ij} + \frac{G}{2\alpha}) +{\rm erfc`}(-\alpha z_{ij} + \frac{G}{2\alpha})\Bigg] \Bigg)\,. \]
Definition at line 1185 of file ion_interaction.F90.
|
private |
Auxiliary function for the Ewald 2D stress.
Definition at line 1328 of file ion_interaction.F90.
subroutine, public ion_interaction_oct_m::ion_interaction_test | ( | class(space_t), intent(in) | space, |
type(lattice_vectors_t), intent(in) | latt, | ||
type(atom_t), dimension(:), intent(in) | atom, | ||
integer, intent(in) | natoms, | ||
real(real64), dimension(1:space%dim,1:natoms), intent(in) | pos, | ||
real(real64), dimension(:), intent(in) | lsize, | ||
type(namespace_t), intent(in) | namespace, | ||
type(multicomm_t), intent(in) | mc | ||
) |
Definition at line 1345 of file ion_interaction.F90.
|
private |
Definition at line 154 of file ion_interaction.F90.
|
private |
Definition at line 154 of file ion_interaction.F90.
|
private |
Definition at line 154 of file ion_interaction.F90.
|
private |
Definition at line 154 of file ion_interaction.F90.