The Elk Code
vecfbz.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2008 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU Lesser General Public
4 ! License. See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: vecfbz
8 ! !INTERFACE:
9 subroutine vecfbz(eps,bvec,vpl)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! eps : zero component tolerance (in,real)
12 ! bvec : reciprocal lattice vectors (in,real(3,3))
13 ! vpl : input vector in lattice coordinates (inout,real(3))
14 ! !DESCRIPTION:
15 ! Maps a vector in lattice coordinates to the first Brillouin zone. This is
16 ! done by first removing its integer components and then adding primitive
17 ! reciprocal lattice vectors until the shortest vector is found.
18 !
19 ! !REVISION HISTORY:
20 ! Created September 2008 (JKD)
21 !EOP
22 !BOC
23 implicit none
24 ! arguments
25 real(8), intent(in) :: eps,bvec(3,3)
26 real(8), intent(inout) :: vpl(3)
27 ! local variables
28 integer i1,i2,i3,j1,j2,j3
29 real(8) v0(3),v1(3),v2(3),v3(3),t1,t2
30 ! map vector to [0,1) interval
31 call r3frac(eps,vpl)
32 v0(1:3)=bvec(1:3,1)*vpl(1)+bvec(1:3,2)*vpl(2)+bvec(1:3,3)*vpl(3)
33 t1=v0(1)**2+v0(2)**2+v0(3)**2
34 j1=0; j2=0; j3=0
35 do i1=-1,0
36  v1(1:3)=v0(1:3)+dble(i1)*bvec(1:3,1)
37  do i2=-1,0
38  v2(1:3)=v1(1:3)+dble(i2)*bvec(1:3,2)
39  do i3=-1,0
40  v3(1:3)=v2(1:3)+dble(i3)*bvec(1:3,3)
41  t2=v3(1)**2+v3(2)**2+v3(3)**2
42  if (t2 < t1+eps) then
43  j1=i1; j2=i2; j3=i3
44  t1=t2
45  end if
46  end do
47  end do
48 end do
49 vpl(1)=vpl(1)+dble(j1)
50 vpl(2)=vpl(2)+dble(j2)
51 vpl(3)=vpl(3)+dble(j3)
52 end subroutine
53 !EOC
54 
subroutine vecfbz(eps, bvec, vpl)
Definition: vecfbz.f90:10
pure subroutine r3frac(eps, v)
Definition: r3frac.f90:10