The Elk Code
 
Loading...
Searching...
No Matches
getevecuv.f90
Go to the documentation of this file.
1
2! Copyright (C) 2019 J. K. Dewhurst and S. Sharma.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine getevecuv(ikp,vpl,u,v)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: ikp
11real(8), intent(in) :: vpl(3)
12complex(8), intent(out) :: u(nstsv,nstsv),v(nstsv,nstsv)
13! local variables
14integer isym,ik
15integer recl,nstsv_
16real(8) vkl_(3),t1
17if (ikp > 0) then
18 ik=ikp
19else
20! find the equivalent k-point number and symmetry which rotates vkl to vpl
21 call findkpt(vpl,isym,ik)
22end if
23! find the record length
24inquire(iolength=recl) vkl_,nstsv_,u,v
25!$OMP CRITICAL(u322)
26open(322,file='EVECUV.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
27read(322,rec=ik) vkl_,nstsv_,u,v
28close(322)
29!$OMP END CRITICAL(u322)
30t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
31if (t1 > epslat) then
32 write(*,*)
33 write(*,'("Error(getevecuv): differing vectors for k-point ",I8)') ik
34 write(*,'(" current : ",3G18.10)') vkl(:,ik)
35 write(*,'(" EVECUV.OUT : ",3G18.10)') vkl_
36 write(*,*)
37 stop
38end if
39if (nstsv /= nstsv_) then
40 write(*,*)
41 write(*,'("Error(getevecuv): differing nstsv for k-point ",I8)') ik
42 write(*,'(" current : ",I8)') nstsv
43 write(*,'(" EVECUV.OUT : ",I8)') nstsv_
44 write(*,*)
45 stop
46end if
47end subroutine
48
subroutine findkpt(vpl, isym, ik)
Definition findkpt.f90:7
subroutine getevecuv(ikp, vpl, u, v)
Definition getevecuv.f90:7
real(8) epslat
Definition modmain.f90:24
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471