Evaluation of the dispersion energy expression
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(structure_type), | intent(in) | :: | mol |
Molecular structure data |
||
real(kind=wp), | intent(in) | :: | trans(:,:) |
Lattice points |
||
real(kind=wp), | intent(in) | :: | cutoff |
Real space cutoff |
||
real(kind=wp), | intent(in) | :: | s9 |
Scaling for dispersion coefficients |
||
real(kind=wp), | intent(in) | :: | a1 |
Scaling parameter for critical radius |
||
real(kind=wp), | intent(in) | :: | a2 |
Offset parameter for critical radius |
||
real(kind=wp), | intent(in) | :: | alp |
Exponent of zero damping function |
||
real(kind=wp), | intent(in) | :: | r4r2(:) |
Expectation values for r4 over r2 operator |
||
real(kind=wp), | intent(in) | :: | c6(:,:) |
C6 coefficients for all atom pairs. |
||
real(kind=wp), | intent(inout) | :: | energy(:) |
Dispersion energy |
subroutine get_atm_dispersion_energy(mol, trans, cutoff, s9, a1, a2, alp, r4r2, & & c6, energy) !> Molecular structure data class(structure_type), intent(in) :: mol !> Lattice points real(wp), intent(in) :: trans(:, :) !> Real space cutoff real(wp), intent(in) :: cutoff !> Scaling for dispersion coefficients real(wp), intent(in) :: s9 !> Scaling parameter for critical radius real(wp), intent(in) :: a1 !> Offset parameter for critical radius real(wp), intent(in) :: a2 !> Exponent of zero damping function real(wp), intent(in) :: alp !> Expectation values for r4 over r2 operator real(wp), intent(in) :: r4r2(:) !> C6 coefficients for all atom pairs. real(wp), intent(in) :: c6(:, :) !> Dispersion energy real(wp), intent(inout) :: energy(:) integer :: iat, jat, kat, izp, jzp, kzp, jtr, ktr real(wp) :: vij(3), vjk(3), vik(3), r2ij, r2jk, r2ik, c6ij, c6jk, c6ik, triple real(wp) :: r0ij, r0jk, r0ik, r0, r1, r2, r3, r5, rr, fdmp, ang real(wp) :: cutoff2, c9, dE cutoff2 = cutoff*cutoff !$omp parallel do schedule(runtime) default(none) reduction(+:energy) & !$omp shared(mol, trans, c6, s9, a1, a2, alp, r4r2, cutoff2) & !$omp private(iat, jat, kat, izp, jzp, kzp, jtr, ktr, vij, vjk, vik, & !$omp& r2ij, r2jk, r2ik, c6ij, c6jk, c6ik, triple, r0ij, r0jk, r0ik, r0, & !$omp& r1, r2, r3, r5, rr, fdmp, ang, c9, dE) do iat = 1, mol%nat izp = mol%id(iat) do jat = 1, iat jzp = mol%id(jat) c6ij = c6(jat, iat) r0ij = a1 * sqrt(3*r4r2(jzp)*r4r2(izp)) + a2 do jtr = 1, size(trans, 2) vij(:) = mol%xyz(:, jat) + trans(:, jtr) - mol%xyz(:, iat) r2ij = vij(1)*vij(1) + vij(2)*vij(2) + vij(3)*vij(3) if (r2ij > cutoff2 .or. r2ij < epsilon(1.0_wp)) cycle do kat = 1, jat kzp = mol%id(kat) c6ik = c6(kat, iat) c6jk = c6(kat, jat) c9 = -s9 * sqrt(abs(c6ij*c6ik*c6jk)) r0ik = a1 * sqrt(3*r4r2(kzp)*r4r2(izp)) + a2 r0jk = a1 * sqrt(3*r4r2(kzp)*r4r2(jzp)) + a2 r0 = r0ij * r0ik * r0jk triple = triple_scale(iat, jat, kat) do ktr = 1, size(trans, 2) vik(:) = mol%xyz(:, kat) + trans(:, ktr) - mol%xyz(:, iat) r2ik = vik(1)*vik(1) + vik(2)*vik(2) + vik(3)*vik(3) if (r2ik > cutoff2 .or. r2ik < epsilon(1.0_wp)) cycle vjk(:) = mol%xyz(:, kat) + trans(:, ktr) - mol%xyz(:, jat) & & - trans(:, jtr) r2jk = vjk(1)*vjk(1) + vjk(2)*vjk(2) + vjk(3)*vjk(3) if (r2jk > cutoff2 .or. r2jk < epsilon(1.0_wp)) cycle r2 = r2ij*r2ik*r2jk r1 = sqrt(r2) r3 = r2 * r1 r5 = r3 * r2 fdmp = 1.0_wp / (1.0_wp + 6.0_wp * (r0 / r1)**(alp / 3.0_wp)) ang = 0.375_wp*(r2ij + r2jk - r2ik)*(r2ij - r2jk + r2ik)& & *(-r2ij + r2jk + r2ik) / r5 + 1.0_wp / r3 rr = ang*fdmp dE = rr * c9 * triple energy(iat) = energy(iat) - dE/3 energy(jat) = energy(jat) - dE/3 energy(kat) = energy(kat) - dE/3 end do end do end do end do end do end subroutine get_atm_dispersion_energy