The Elk Code
Loading...
Searching...
No Matches
findkpt.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3
! This file is distributed under the terms of the GNU General Public License.
4
! See the file COPYING for license details.
5
6
subroutine
findkpt
(vpl,isym,ik)
7
use
modmain
8
implicit none
9
! arguments
10
real
(8),
intent(in)
:: vpl(3)
11
integer
,
intent(out)
:: isym,ik
12
! local variables
13
integer
i1,i2,i3,lspl
14
real
(8) v1(3),v2(3),t1
15
v1(1:3)=vpl(1:3)-
vkloff
(1:3)/dble(
ngridk
(1:3))
16
i1=modulo(nint(v1(1)*
ngridk
(1)),
ngridk
(1))
17
i2=modulo(nint(v1(2)*
ngridk
(2)),
ngridk
(2))
18
i3=modulo(nint(v1(3)*
ngridk
(3)),
ngridk
(3))
19
ik=
ivkik
(i1,i2,i3)
20
v1(1:3)=
vkl
(1:3,ik)
21
! find the symmetry which rotates vkl to vpl
22
do
isym=1,
nsymcrys
23
lspl=
lsplsymc
(isym)
24
! multiply vpl by the transpose of the symmetry matrix (i.e. the inverse)
25
v2(1:3)=
symlat
(1,1:3,lspl)*vpl(1) &
26
+
symlat
(2,1:3,lspl)*vpl(2) &
27
+
symlat
(3,1:3,lspl)*vpl(3)
28
call
r3frac
(
epslat
,v2)
29
t1=abs(v1(1)-v2(1))+abs(v1(2)-v2(2))+abs(v1(3)-v2(3))
30
if
(t1 <
epslat
)
return
31
end do
32
write
(*,*)
33
write
(*,
'("Error(findkpt): equivalent k-point not in set")'
)
34
write
(*,
'(" Requested k-point : ",3G18.10)'
) vpl
35
write
(*,*)
36
stop
37
end subroutine
38
findkpt
subroutine findkpt(vpl, isym, ik)
Definition
findkpt.f90:7
modmain
Definition
modmain.f90:6
modmain::epslat
real(8) epslat
Definition
modmain.f90:24
modmain::ngridk
integer, dimension(3) ngridk
Definition
modmain.f90:448
modmain::vkl
real(8), dimension(:,:), allocatable vkl
Definition
modmain.f90:471
modmain::ivkik
integer, dimension(:,:,:), allocatable ivkik
Definition
modmain.f90:467
modmain::vkloff
real(8), dimension(3) vkloff
Definition
modmain.f90:450
modmain::nsymcrys
integer nsymcrys
Definition
modmain.f90:358
modmain::symlat
integer, dimension(3, 3, 48) symlat
Definition
modmain.f90:344
modmain::lsplsymc
integer, dimension(maxsymcrys) lsplsymc
Definition
modmain.f90:364
r3frac
pure subroutine r3frac(eps, v)
Definition
r3frac.f90:10
findkpt.f90
Generated by
1.9.8