The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine getdevecsv(ik,iq,is,ia,ip,devecsv)
7use modmain
8use modphonon
10implicit none
11! arguments
12integer, intent(in) :: ik,iq,is,ia,ip
13complex(8), intent(out) :: devecsv(nstsv,nstsv)
14! local variables
15logical tgs
16integer recl,nstsv_
17real(8) vkl_(3),t1
18character(256) fext,fname
19! construct the dynamical matrix file extension
20call dynfext(iq,is,ia,ip,fext)
21! construct filename
22fname=trim(scrpath)//'DEVECSV'//trim(fext)
23!$OMP CRITICAL(u226)
24! read from RAM disk if required
25if (ramdisk) then
26 call getrd(fname,ik,tgs,v1=vkl_,n1=nstsv_,nzv=nstsv*nstsv,zva=devecsv)
27 if (tgs) goto 10
28end if
29! find the record length
30inquire(iolength=recl) vkl_,nstsv_,devecsv
31open(226,file=fname,form='UNFORMATTED',access='DIRECT',recl=recl)
32read(226,rec=ik) vkl_,nstsv_,devecsv
33close(226)
3410 continue
35!$OMP END CRITICAL(u226)
36t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
37if (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
44end if
45if (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
52end if
53end subroutine
54
subroutine dynfext(iq, is, ia, ip, fext)
Definition dynfext.f90:7
subroutine getdevecsv(ik, iq, is, ia, ip, devecsv)
Definition getdevecsv.f90:7
character(256) scrpath
Definition modmain.f90:1302
real(8) epslat
Definition modmain.f90:24
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
type(file_t), dimension(:), allocatable, private file
subroutine getrd(fname, irec, tgs, n1, n2, n3, v1, v2, nrv, rva, nzv, zva)
logical ramdisk
Definition modramdisk.f90:9