The Elk Code
rfhkintp.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2015 D. Ernsting, S. Dugdale and J. K. Dewhurst.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 real(8) function rfhkintp(vhpl,rfhk)
7 use modmain
8 use modpw
9 implicit none
10 ! arguments
11 real(8), intent(in) :: vhpl(3)
12 real(4), intent(in) :: rfhk(nhkmax,nkpt)
13 ! local variables
14 integer ivh0(3),ivk0(3),ihk
15 integer ivhb(3,0:1,0:1,0:1)
16 integer ivkb(3,0:1,0:1,0:1)
17 integer isym,lspl,ik,jk,i,j,k
18 real(8) vpl(3),fb(0:1,0:1,0:1)
19 real(8) f00,f01,f10,f11,f0,f1
20 real(8) v0(3),v1(3),v2(3),t1,t2
21 ! find the H-vector and k-vector corresponding to the input H+p-vector
22 ivh0(:)=floor(vhpl(:))
23 vpl(:)=vhpl(:)-dble(ivh0(:))
24 v1(:)=vpl(:)*dble(ngridk(:))
25 ivk0(:)=floor(v1(:))
26 ! determine the corners of the box containing the input point
27 do i=0,1; do j=0,1; do k=0,1
28  ivkb(1,i,j,k)=ivk0(1)+i
29  ivkb(2,i,j,k)=ivk0(2)+j
30  ivkb(3,i,j,k)=ivk0(3)+k
31  ivhb(:,i,j,k)=ivh0(:)+ivkb(:,i,j,k)/ngridk(:)
32  ivkb(:,i,j,k)=modulo(ivkb(:,i,j,k),ngridk(:))
33 end do; end do; end do
34 ! determine the function at each corner of the box
35 do i=0,1; do j=0,1; do k=0,1
36  fb(i,j,k)=0.d0
37 ! non-reduced k-point index
38  jk=ivkiknr(ivkb(1,i,j,k),ivkb(2,i,j,k),ivkb(3,i,j,k))
39 ! H+k-vector at corner of box
40  v1(:)=dble(ivhb(:,i,j,k))+vkl(:,jk)
41 ! store the origin of the box
42  if ((i == 0).and.(j == 0).and.(k == 0)) v0(:)=v1(:)
43 ! vector in Cartesian coordinates
44  v2(:)=bvec(:,1)*v1(1)+bvec(:,2)*v1(2)+bvec(:,3)*v1(3)
45 ! check length is within range
46  t1=sqrt(v2(1)**2+v2(2)**2+v2(3)**2)
47  if (t1 > hkmax) cycle
48 ! find the lattice symmetry which maps the non-reduced to reduced k-point
49  call findkpt(vkl(:,jk),isym,ik)
50 ! index to spatial rotation in lattice point group
51  lspl=lsplsymc(isym)
52  v2(:)=symlat(1,:,lspl)*v1(1)+symlat(2,:,lspl)*v1(2)+symlat(3,:,lspl)*v1(3)
53 ! find the H+k-vector for the reduced k-point
54  do ihk=1,nhk(1,ik)
55  t1=abs(v2(1)-vhkl(1,ihk,1,ik)) &
56  +abs(v2(2)-vhkl(2,ihk,1,ik)) &
57  +abs(v2(3)-vhkl(3,ihk,1,ik))
58  if (t1 < epslat) then
59  fb(i,j,k)=rfhk(ihk,ik)
60  exit
61  end if
62  end do
63 end do; end do; end do
64 ! interpolate function
65 t2=(vhpl(1)-v0(1))*dble(ngridk(1))
66 t1=1.d0-t2
67 f00=fb(0,0,0)*t1+fb(1,0,0)*t2
68 f01=fb(0,0,1)*t1+fb(1,0,1)*t2
69 f10=fb(0,1,0)*t1+fb(1,1,0)*t2
70 f11=fb(0,1,1)*t1+fb(1,1,1)*t2
71 t2=(vhpl(2)-v0(2))*dble(ngridk(2))
72 t1=1.d0-t2
73 f0=f00*t1+f10*t2
74 f1=f01*t1+f11*t2
75 t2=(vhpl(3)-v0(3))*dble(ngridk(3))
76 t1=1.d0-t2
77 rfhkintp=f0*t1+f1*t2
78 end function
79 
integer nkpt
Definition: modmain.f90:461
integer, dimension(3, 3, 48) symlat
Definition: modmain.f90:344
integer, dimension(maxsymcrys) lsplsymc
Definition: modmain.f90:364
integer, dimension(3) ngridk
Definition: modmain.f90:448
real(8), dimension(3, 3) bvec
Definition: modmain.f90:16
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
real(8) epslat
Definition: modmain.f90:24
Definition: modpw.f90:6
real(8) function rfhkintp(vhpl, rfhk)
Definition: rfhkintp.f90:7
subroutine findkpt(vpl, isym, ik)
Definition: findkpt.f90:7
integer nhkmax
Definition: modpw.f90:42
integer, dimension(:,:,:), allocatable ivkiknr
Definition: modmain.f90:469