The Elk Code
 
Loading...
Searching...
No Matches
getdevecfv.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 getdevecfv(ik,iq,is,ia,ip,devecfv)
7use modmain
8use modphonon
10implicit none
11! arguments
12integer, intent(in) :: ik,iq,is,ia,ip
13complex(8), intent(out) :: devecfv(nmatmax,nstfv,nspnfv)
14! local variables
15logical tgs
16integer recl,nmatmax_,nstfv_,nspnfv_
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)//'DEVECFV'//trim(fext)
23!$OMP CRITICAL(u222)
24! read from RAM disk if required
25if (ramdisk) then
26 call getrd(fname,ik,tgs,v1=vkl_,n1=nmatmax_,n2=nstfv_,n3=nspnfv_, &
27 nzv=nmatmax*nstfv*nspnfv,zva=devecfv)
28 if (tgs) goto 10
29end if
30! find the record length
31inquire(iolength=recl) vkl_,nmatmax_,nstfv_,nspnfv_,devecfv
32open(222,file=fname,form='UNFORMATTED',access='DIRECT',recl=recl)
33read(222,rec=ik) vkl_,nmatmax_,nstfv_,nspnfv_,devecfv
34close(222)
3510 continue
36!$OMP END CRITICAL(u222)
37t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
38if (t1 > epslat) then
39 write(*,*)
40 write(*,'("Error(getdevecfv): differing vectors for k-point ",I8)') ik
41 write(*,'(" current : ",3G18.10)') vkl(:,ik)
42 write(*,'(" ",A," : ",3G18.10)') trim(fname),vkl_
43 write(*,*)
44 stop
45end if
46if (nmatmax /= nmatmax_) then
47 write(*,*)
48 write(*,'("Error(getdevecfv): differing nmatmax for k-point ",I8)') ik
49 write(*,'(" current : ",I8)') nmatmax
50 write(*,'(" ",A," : ",I8)') trim(fname),nmatmax_
51 write(*,*)
52 stop
53end if
54if (nstfv /= nstfv_) then
55 write(*,*)
56 write(*,'("Error(getdevecfv): differing nstfv for k-point ",I8)') ik
57 write(*,'(" current : ",I8)') nstfv
58 write(*,'(" ",A," : ",I8)') trim(fname),nstfv_
59 write(*,*)
60 stop
61end if
62if (nspnfv /= nspnfv_) then
63 write(*,*)
64 write(*,'("Error(getdevecfv): differing nspnfv for k-point ",I8)') ik
65 write(*,'(" current : ",I8)') nspnfv
66 write(*,'(" ",A," : ",I8)') trim(fname),nspnfv_
67 write(*,*)
68 stop
69end if
70end subroutine
71
subroutine dynfext(iq, is, ia, ip, fext)
Definition dynfext.f90:7
subroutine getdevecfv(ik, iq, is, ia, ip, devecfv)
Definition getdevecfv.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