The Elk Code
getevalfv.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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 getevalfv(fext,ikp,vpl,evalfv)
7 use modmain
8 use modramdisk
9 implicit none
10 ! arguments
11 character(*), intent(in) :: fext
12 integer, intent(in) :: ikp
13 real(8), intent(in) :: vpl(3)
14 real(8), intent(out) :: evalfv(nstfv,nspnfv)
15 ! local variables
16 logical tgs
17 integer isym,ik
18 integer recl,nstfv_,nspnfv_
19 real(8) vkl_(3),t1
20 character(256) fname
21 if (ikp > 0) then
22  ik=ikp
23 else
24 ! find the k-point number
25  call findkpt(vpl,isym,ik)
26 end if
27 ! construct the filename
28 fname='EVALFV'//trim(fext)
29 !$OMP CRITICAL(u200)
30 ! read from RAM disk if required
31 if (ramdisk) then
32  call getrd(fname,ik,tgs,v1=vkl_,n1=nstfv_,n2=nspnfv_,nrv=nstfv*nspnfv, &
33  rva=evalfv)
34  if (tgs) goto 10
35 end if
36 ! find the record length
37 inquire(iolength=recl) vkl_,nstfv_,nspnfv_,evalfv
38 open(200,file=fname,form='UNFORMATTED',access='DIRECT',recl=recl)
39 read(200,rec=ik) vkl_,nstfv_,nspnfv_,evalfv
40 close(200)
41 10 continue
42 !$OMP END CRITICAL(u200)
43 t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
44 if (t1 > epslat) then
45  write(*,*)
46  write(*,'("Error(getevalfv): differing vectors for k-point ",I8)') ik
47  write(*,'(" current : ",3G18.10)') vkl(:,ik)
48  write(*,'(" EVALFV.OUT : ",3G18.10)') vkl_
49  write(*,*)
50  stop
51 end if
52 if (nstfv /= nstfv_) then
53  write(*,*)
54  write(*,'("Error(getevalfv): differing nstfv for k-point ",I8)') ik
55  write(*,'(" current : ",I8)') nstfv
56  write(*,'(" EVALFV.OUT : ",I8)') nstfv_
57  write(*,*)
58  stop
59 end if
60 if (nspnfv /= nspnfv_) then
61  write(*,*)
62  write(*,'("Error(getevalfv): differing nspnfv for k-point ",I8)') ik
63  write(*,'(" current : ",I8)') nspnfv
64  write(*,'(" EVALFV.OUT : ",I8)') nspnfv_
65  write(*,*)
66  stop
67 end if
68 end subroutine
69 
logical ramdisk
Definition: modramdisk.f90:9
type(file_t), dimension(:), allocatable, private file
Definition: modramdisk.f90:29
subroutine getrd(fname, irec, tgs, n1, n2, n3, v1, v2, nrv, rva, nzv, zva)
Definition: modramdisk.f90:214
subroutine getevalfv(fext, ikp, vpl, evalfv)
Definition: getevalfv.f90:7
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
real(8) epslat
Definition: modmain.f90:24
subroutine findkpt(vpl, isym, ik)
Definition: findkpt.f90:7