The Elk Code
getdevecsv.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2013 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 
6 subroutine getdevecsv(ik,iq,is,ia,ip,devecsv)
7 use modmain
8 use modphonon
9 use modramdisk
10 implicit none
11 ! arguments
12 integer, intent(in) :: ik,iq,is,ia,ip
13 complex(8), intent(out) :: devecsv(nstsv,nstsv)
14 ! local variables
15 logical tgs
16 integer recl,nstsv_
17 real(8) vkl_(3),t1
18 character(256) fext,fname
19 ! construct the dynamical matrix file extension
20 call dynfext(iq,is,ia,ip,fext)
21 ! construct filename
22 fname=trim(scrpath)//'DEVECSV'//trim(fext)
23 !$OMP CRITICAL(u226)
24 ! read from RAM disk if required
25 if (ramdisk) then
26  call getrd(fname,ik,tgs,v1=vkl_,n1=nstsv_,nzv=nstsv*nstsv,zva=devecsv)
27  if (tgs) goto 10
28 end if
29 ! find the record length
30 inquire(iolength=recl) vkl_,nstsv_,devecsv
31 open(226,file=fname,form='UNFORMATTED',access='DIRECT',recl=recl)
32 read(226,rec=ik) vkl_,nstsv_,devecsv
33 close(226)
34 10 continue
35 !$OMP END CRITICAL(u226)
36 t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
37 if (t1 > epslat) then
38  write(*,*)
39  write(*,'("Error(getdevecsv): differing vectors for k-point ",I8)') ik
40  write(*,'(" current : ",3G18.10)') vkl(:,ik)
41  write(*,'(" ",A," : ",3G18.10)') trim(fname),vkl_
42  write(*,*)
43  stop
44 end if
45 if (nstsv /= nstsv_) then
46  write(*,*)
47  write(*,'("Error(getdevecsv): differing nstsv for k-point ",I8)') ik
48  write(*,'(" current : ",I8)') nstsv
49  write(*,'(" ",A," : ",I8)') trim(fname),nstsv_
50  write(*,*)
51  stop
52 end if
53 end subroutine
54 
character(256) scrpath
Definition: modmain.f90:1303
logical ramdisk
Definition: modramdisk.f90:9
subroutine dynfext(iq, is, ia, ip, fext)
Definition: dynfext.f90:7
type(file_t), dimension(:), allocatable, private file
Definition: modramdisk.f90:29
subroutine getdevecsv(ik, iq, is, ia, ip, devecsv)
Definition: getdevecsv.f90:7
subroutine getrd(fname, irec, tgs, n1, n2, n3, v1, v2, nrv, rva, nzv, zva)
Definition: modramdisk.f90:214
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
real(8) epslat
Definition: modmain.f90:24