The Elk Code
 
Loading...
Searching...
No Matches
findqpt.f90
Go to the documentation of this file.
1
2! Copyright (C) 2010 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
6subroutine findqpt(vpl,isym,iq)
7use modmain
8implicit none
9! arguments
10real(8), intent(in) :: vpl(3)
11integer, intent(out) :: isym,iq
12! local variables
13integer i1,i2,i3,lspl
14real(8) v1(3),v2(3),t1
15i1=modulo(nint(vpl(1)*ngridq(1)),ngridq(1))
16i2=modulo(nint(vpl(2)*ngridq(2)),ngridq(2))
17i3=modulo(nint(vpl(3)*ngridq(3)),ngridq(3))
18iq=ivqiq(i1,i2,i3)
19v1(1:3)=vql(1:3,iq)
20call r3frac(epslat,v1)
21! find the symmetry which rotates vql to vpl
22do 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
31end do
32write(*,*)
33write(*,'("Error(findqpt): equivalent q-point not in set")')
34write(*,'(" Requested q-point : ",3G18.10)') vpl
35write(*,*)
36stop
37end subroutine
38
subroutine findqpt(vpl, isym, iq)
Definition findqpt.f90:7
integer, dimension(:,:,:), allocatable ivqiq
Definition modmain.f90:531
real(8) epslat
Definition modmain.f90:24
real(8), dimension(:,:), allocatable vql
Definition modmain.f90:545
integer nsymcrys
Definition modmain.f90:358
integer, dimension(3, 3, 48) symlat
Definition modmain.f90:344
integer, dimension(3) ngridq
Definition modmain.f90:515
integer, dimension(maxsymcrys) lsplsymc
Definition modmain.f90:364
pure subroutine r3frac(eps, v)
Definition r3frac.f90:10