Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
96 changes: 63 additions & 33 deletions src/dftd3/cutoff.f90
Original file line number Diff line number Diff line change
Expand Up @@ -55,69 +55,99 @@ subroutine get_lattice_points(periodic, lat, rthr, trans)
real(wp), intent(in) :: rthr
real(wp), intent(in) :: lat(:, :)
real(wp), allocatable, intent(out) :: trans(:, :)
real(wp) :: vec(size(lat, 2)), norms(size(lat, 1), size(lat, 2))
integer :: rep(3)
integer :: itr, ix, iy, iz

if (.not.any(periodic)) then
allocate(trans(3, 1))
trans(:, :) = 0.0_wp
else
return
end if

if (all(periodic)) then
call get_translations(lat, rthr, rep)
allocate(trans(3, product(2*rep+1)))
itr = 0
do ix = -rep(1), rep(1)
do iy = -rep(2), rep(2)
do iz = -rep(3), rep(3)
itr = itr + 1
trans(:, itr) = lat(:, 1)*ix + lat(:, 2)*iy + lat(:, 3)*iz
end do
else if (count(periodic) == 1) then
vec = sqrt(sum(lat**2, 1))
where(periodic)
rep = ceiling(rthr / vec)
elsewhere
rep = 0
endwhere
else if (count(periodic) == 2) then
call get_normals(lat, norms)
where(spread(periodic, 2, 3))
norms = lat
endwhere
call get_translations(norms, rthr, rep)
where(.not.periodic)
rep = 0
endwhere
end if

allocate(trans(3, product(2*rep+1)))
itr = 0
do ix = -rep(1), rep(1)
do iy = -rep(2), rep(2)
do iz = -rep(3), rep(3)
itr = itr + 1
trans(:, itr) = lat(:, 1)*ix + lat(:, 2)*iy + lat(:, 3)*iz
end do
end do
end if
end do

end subroutine get_lattice_points


subroutine get_normals(lattice, normal)
real(wp), intent(in) :: lattice(:, :)
real(wp), intent(out) :: normal(:, :)

integer :: itr
do itr = 1, 3
call crossproduct(lattice(:, mod(itr, 3) + 1), lattice(:, mod(itr + 1, 3) + 1), &
& normal(:, itr))
end do
end subroutine get_normals


!> generate a supercell based on a realspace cutoff, this subroutine
! doesn't know anything about the convergence behaviour of the
! associated property.
pure subroutine get_translations(lat, rthr, rep)
real(wp), intent(in) :: rthr
real(wp), intent(in) :: lat(3, 3)
integer, intent(out) :: rep(3)
real(wp) :: normx(3), normy(3), normz(3)
real(wp) :: norm(3, 3), normy(3), normz(3)
real(wp) :: cos10, cos21, cos32

! find normal to the plane...
call crossproduct(lat(:, 2), lat(:, 3), normx)
call crossproduct(lat(:, 3), lat(:, 1), normy)
call crossproduct(lat(:, 1), lat(:, 2), normz)
call crossproduct(lat(:, 2), lat(:, 3), norm(:, 1))
call crossproduct(lat(:, 3), lat(:, 1), norm(:, 2))
call crossproduct(lat(:, 1), lat(:, 2), norm(:, 3))
! ...normalize it...
normx = normx/norm2(normx)
normy = normy/norm2(normy)
normz = normz/norm2(normz)
norm(:, 1) = norm(:, 1)/norm2(norm(:, 1))
norm(:, 2) = norm(:, 2)/norm2(norm(:, 2))
norm(:, 3) = norm(:, 3)/norm2(norm(:, 3))
! cos angles between normals and lattice vectors
cos10 = sum(normx*lat(:, 1))
cos21 = sum(normy*lat(:, 2))
cos32 = sum(normz*lat(:, 3))
cos10 = sum(norm(:, 1)*lat(:, 1))
cos21 = sum(norm(:, 2)*lat(:, 2))
cos32 = sum(norm(:, 3)*lat(:, 3))
rep(1) = ceiling(abs(rthr/cos10))
rep(2) = ceiling(abs(rthr/cos21))
rep(3) = ceiling(abs(rthr/cos32))

contains

pure subroutine crossproduct(a, b, c)
real(wp), intent(in) :: a(3), b(3)
real(wp), intent(out) :: c(3)
real(wp) :: x, y, z

x=a(2)*b(3)-b(2)*a(3)
y=a(3)*b(1)-b(3)*a(1)
z=a(1)*b(2)-b(1)*a(2)
c=(/x, y, z/)
end subroutine crossproduct

end subroutine get_translations

pure subroutine crossproduct(a, b, c)
real(wp), intent(in) :: a(3), b(3)
real(wp), intent(out) :: c(3)
real(wp) :: x, y, z

x = a(2)*b(3)-b(2)*a(3)
y = a(3)*b(1)-b(3)*a(1)
z = a(1)*b(2)-b(1)*a(2)
c = [x, y, z]
end subroutine crossproduct

end module dftd3_cutoff
108 changes: 106 additions & 2 deletions test/unit/test_periodic.f90
Original file line number Diff line number Diff line change
Expand Up @@ -18,7 +18,7 @@ module test_periodic
use mctc_env, only : wp
use mctc_env_testing, only : new_unittest, unittest_type, error_type, check, &
& test_failed
use mctc_io, only : structure_type
use mctc_io, only : structure_type, new
use mstore, only : get_structure
use dftd3
implicit none
Expand Down Expand Up @@ -50,11 +50,13 @@ subroutine collect_periodic(testsuite)
& new_unittest("BLYP-D3(0)", test_blypd3zero_benzene), &
& new_unittest("M06L-D3(0)", test_m06ld3zero_cyanamide), &
& new_unittest("rPW86PBE-D3(0)", test_rpw86pbed3zero_co2), &
& new_unittest("revSSB-D3(0)", test_revssbd3zero_cytosine) &
& new_unittest("revSSB-D3(0)", test_revssbd3zero_cytosine), &
! new_unittest("HSEsol-D3(BJ)-ATM", test_hsesold3bjatm_oxacb), &
! new_unittest("PWGGA-D3(BJ)-ATM", test_pwggad3bjatm_pyrazine), &
! new_unittest("B3PW91-D3(0)-ATM", test_b3pw91d3zeroatm_urea), &
! new_unittest("RPBE-D3(0)-ATM", test_rpbed3zeroatm_hexamine) &
& new_unittest("revPBE-D3(BJ)", test_revpbed3bj_1d), &
& new_unittest("BP86-D3(0)", test_bp86d3zero_2d) &
& ]

end subroutine collect_periodic
Expand Down Expand Up @@ -397,4 +399,106 @@ subroutine test_rpbed3zeroatm_hexamine(error)
end subroutine test_rpbed3zeroatm_hexamine


subroutine test_revpbed3bj_1d(error)
!> Error handling
type(error_type), allocatable, intent(out) :: error

integer, parameter :: nat = 32
integer, parameter :: num(nat) = 6
real(wp), parameter :: xyz(3, nat) = reshape([&
& 1.36794785746435_wp, 13.45808943446053_wp, 8.83754983226359_wp, &
& 3.69183290816438_wp, 13.13552229161569_wp, 10.16652201690950_wp, &
& 1.36792668081267_wp, 10.38660504434782_wp, 13.04411926632965_wp, &
& 3.69180534781206_wp, 11.55414582295511_wp, 12.33193380846742_wp, &
& 1.36791549262702_wp, 3.53066844289674_wp, 10.38660588677206_wp, &
& 1.36792046664920_wp, 7.73723910626293_wp, 13.45809224934817_wp, &
& 3.69181279359489_wp, 6.40826717723392_wp, 13.13552570942280_wp, &
& 1.36792009865062_wp, 3.11669712338516_wp, 7.73723850632628_wp, &
& 3.69181515738094_wp, 3.43926499914873_wp, 6.40826580885474_wp, &
& 3.69178443989294_wp, 4.24285720771059_wp, 11.55415026712869_wp, &
& 1.36790824853106_wp, 6.18818490375705_wp, 3.53066863732142_wp, &
& 3.69178194163078_wp, 5.02063901427657_wp, 4.24285736953327_wp, &
& 1.36794124909207_wp, 13.04411858182861_wp, 6.18818324080182_wp, &
& 1.36792249732236_wp, 8.83755133592807_wp, 3.11669686076913_wp, &
& 3.69182456413952_wp, 10.16652118921143_wp, 3.43926084011816_wp, &
& 3.69181444966104_wp, 12.33193631088573_wp, 5.02063847821044_wp, &
& 6.01572566324028_wp, 13.45790756713123_wp, 8.83752222635545_wp, &
& 8.33965926123256_wp, 13.13576644753615_wp, 10.16660228658307_wp, &
& 6.01574747573805_wp, 10.38654070512969_wp, 13.04391961251944_wp, &
& 8.33964066450677_wp, 11.55427002850905_wp, 12.33211653730939_wp, &
& 6.01574728097580_wp, 3.53087013230607_wp, 10.38654217813321_wp, &
& 6.01568913853645_wp, 7.73726406411719_wp, 13.45790864082374_wp, &
& 8.33963586549168_wp, 6.40818371470975_wp, 13.13576911116618_wp, &
& 6.01568179676984_wp, 3.11688332536281_wp, 7.73726611148835_wp, &
& 8.33963704688671_wp, 3.43902559351770_wp, 6.40818390180453_wp, &
& 8.33962496288127_wp, 4.24267007149867_wp, 11.55427031066552_wp, &
& 6.01573464280675_wp, 6.18824653544318_wp, 3.53086861480278_wp, &
& 8.33961857277245_wp, 5.02052001792996_wp, 4.24267413625204_wp, &
& 6.01575677304189_wp, 13.04392044501564_wp, 6.18824448603611_wp, &
& 6.01568344836224_wp, 8.83752193432504_wp, 3.11688171781516_wp, &
& 8.33964228963694_wp, 10.16660428027860_wp, 3.43902155668011_wp, &
& 8.33965118613331_wp, 12.33211762632282_wp, 5.02051902430387_wp],&
& shape(xyz))
logical, parameter :: periodic(3) = [.true., .false., .false.]
real(wp), parameter :: lattice(3, 3) = reshape([&
& 9.2955628527586_wp, 0.0_wp, 0.0_wp, 0.0_wp, 1.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 1.0_wp], &
& shape(lattice))
real(wp), parameter :: lattice_vac(3, 3) = reshape([&
& 9.2955628527586_wp, 0.0_wp, 0.0_wp, 0.0_wp, 100.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 100.0_wp], &
& shape(lattice))
type(structure_type) :: mol
type(rational_damping_param) :: param
type(d3_param) :: inp = d3_param(&
& s6 = 1.0_wp, s9 = 0.0_wp, alp = 14.0_wp, &
& a1 = 0.5238_wp, s8 = 2.3550_wp, a2 = 3.5016_wp)

call new(mol, num, xyz, periodic=periodic, lattice=lattice)
call new_rational_damping(param, inp)
call test_dftd3_gen(error, mol, param, -0.25436101196186600_wp)

call new(mol, num, xyz, periodic=[.true.], lattice=lattice_vac)
call test_dftd3_gen(error, mol, param, -0.25436101196186600_wp)

end subroutine test_revpbed3bj_1d


subroutine test_bp86d3zero_2d(error)
!> Error handling
type(error_type), allocatable, intent(out) :: error

integer, parameter :: nat = 4
integer, parameter :: num(nat) = 6
real(wp), parameter :: xyz(3, nat) = reshape([&
& -0.12918412100093_wp, 0.06210659750976_wp, -2.13384498734326_wp, &
& 0.12856915667443_wp, -0.07403227791901_wp, 4.02358027265954_wp, &
& -0.12317720857511_wp, 2.75170732207802_wp, -2.13345350602279_wp, &
& 2.44816466162280_wp, 1.28612566399214_wp, 4.02317048854901_wp],&
& shape(xyz))
logical, parameter :: periodic(3) = [.true., .true., .false.]
real(wp), parameter :: lattice(3, 3) = reshape([&
& 4.68837849314507_wp, 0.00000000000000_wp, 0.00000000000000_wp, &
& -2.36282788044783_wp, 4.04978545156612_wp, 0.00000000000000_wp, &
& 0.00000000000000_wp, 0.00000000000000_wp, 1.00000000000000_wp],&
& shape(lattice))
real(wp), parameter :: lattice_vac(3, 3) = reshape([&
& 4.68837849314507_wp, 0.00000000000000_wp, 0.00000000000000_wp, &
& -2.36282788044783_wp, 4.04978545156612_wp, 0.00000000000000_wp, &
& 0.00000000000000_wp, 0.00000000000000_wp, 100.000000000000_wp],&
& shape(lattice))
type(structure_type) :: mol
type(zero_damping_param) :: param
type(d3_param) :: inp = d3_param(&
& s6 = 1.0_wp, s9 = 0.0_wp, alp = 14.0_wp, &
& rs6 = 1.139_wp, s8 = 1.683_wp)

call new(mol, num, xyz, periodic=periodic, lattice=lattice)
call new_zero_damping(param, inp)
call test_dftd3_gen(error, mol, param, -1.5983949078995030E-2_wp)

call new(mol, num, xyz, periodic=[.true.], lattice=lattice_vac)
call test_dftd3_gen(error, mol, param, -1.5983949078995030E-2_wp)

end subroutine test_bp86d3zero_2d


end module test_periodic