The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine getevalfv(fext,ikp,vpl,evalfv)
7use modmain
9implicit none
10! arguments
11character(*), intent(in) :: fext
12integer, intent(in) :: ikp
13real(8), intent(in) :: vpl(3)
14real(8), intent(out) :: evalfv(nstfv,nspnfv)
15! local variables
16logical tgs
17integer isym,ik
18integer recl,nstfv_,nspnfv_
19real(8) vkl_(3),t1
20character(256) fname
21if (ikp > 0) then
22 ik=ikp
23else
24! find the k-point number
25 call findkpt(vpl,isym,ik)
26end if
27! construct the filename
28fname='EVALFV'//trim(fext)
29!$OMP CRITICAL(u200)
30! read from RAM disk if required
31if (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
35end if
36! find the record length
37inquire(iolength=recl) vkl_,nstfv_,nspnfv_,evalfv
38open(200,file=fname,form='UNFORMATTED',access='DIRECT',recl=recl)
39read(200,rec=ik) vkl_,nstfv_,nspnfv_,evalfv
40close(200)
4110 continue
42!$OMP END CRITICAL(u200)
43t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
44if (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
51end if
52if (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
59end if
60if (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
67end if
68end subroutine
69
subroutine findkpt(vpl, isym, ik)
Definition findkpt.f90:7
subroutine getevalfv(fext, ikp, vpl, evalfv)
Definition getevalfv.f90:7
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