The Elk Code
mtdmin.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2011 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3
! This file is distributed under the terms of the GNU General Public License.
4
! See the file COPYING for license details.
5
6
!BOP
7
! !ROUTINE: mtdmin
8
! !INTERFACE:
9
pure subroutine
mtdmin
(is,js,dmin)
10
! !USES:
11
use
modmain
12
! !INPUT/OUTPUT PARAMETERS:
13
! is, js : species numbers (out,integer)
14
! dmin : minimum distance between muffin-tin surfaces (out,real)
15
! !DESCRIPTION:
16
! Finds the atomic species pair for which the distance between the muffin-tin
17
! surfaces is a minimum. This distance may be negative if the muffin-tins
18
! overlap.
19
!
20
! !REVISION HISTORY:
21
! Created October 2011 (JKD)
22
!EOP
23
!BOC
24
implicit none
25
! arguments
26
integer
,
intent(out)
:: is,js
27
real(8)
,
intent(out)
:: dmin
28
! local variables
29
integer
i1,i2,i3,ks,ka,ls,la
30
real(8)
d2,d2m,r,r2,rm,dr,drm
31
real(8)
v1,v2,v3,w1,w2,w3
32
is=1
33
js=1
34
drm=1.d8
35
do
i1=-1,1;
do
i2=-1,1;
do
i3=-1,1
36
v1=i1*
avec
(1,1)+i2*
avec
(1,2)+i3*
avec
(1,3)
37
v2=i1*
avec
(2,1)+i2*
avec
(2,2)+i3*
avec
(2,3)
38
v3=i1*
avec
(3,1)+i2*
avec
(3,2)+i3*
avec
(3,3)
39
do
ks=1,
nspecies
40
do
ka=1,
natoms
(ks)
41
w1=v1+
atposc
(1,ka,ks)
42
w2=v2+
atposc
(2,ka,ks)
43
w3=v3+
atposc
(3,ka,ks)
44
do
ls=1,
nspecies
45
r=
rmt
(ks)+
rmt
(ls)
46
r2=r**2
47
do
la=1,
natoms
(ls)
48
if
((i1 /= 0).or.(i2 /= 0).or.(i3 /= 0).or. &
49
(ks /= ls).or.(ka /= la))
then
50
d2=(w1-
atposc
(1,la,ls))**2 &
51
+(w2-
atposc
(2,la,ls))**2 &
52
+(w3-
atposc
(3,la,ls))**2
53
dr=d2-r2
54
if
(dr < drm-
epslat
)
then
55
is=ks
56
js=ls
57
rm=r
58
d2m=d2
59
drm=dr
60
end if
61
end if
62
end do
63
end do
64
end do
65
end do
66
end
do; end do; end do
67
dmin=sqrt(d2m)-rm
68
end subroutine
69
!EOC
70
mtdmin
pure subroutine mtdmin(is, js, dmin)
Definition:
mtdmin.f90:10
modmain
Definition:
modmain.f90:6
modmain::avec
real(8), dimension(3, 3) avec
Definition:
modmain.f90:12
modmain::rmt
real(8), dimension(maxspecies) rmt
Definition:
modmain.f90:162
modmain::natoms
integer, dimension(maxspecies) natoms
Definition:
modmain.f90:36
modmain::epslat
real(8) epslat
Definition:
modmain.f90:24
modmain::nspecies
integer nspecies
Definition:
modmain.f90:34
modmain::atposc
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition:
modmain.f90:54
mtdmin.f90
Generated by
1.8.14