The Elk Code
getephmat.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2019 Chung-Yu Wang, J. K. Dewhurst, S. Sharma and
3 ! E. K. U. Gross. This file is distributed under the terms of the GNU General
4 ! Public License. See the file COPYING for license details.
5 
6 subroutine getephmat(iqp,ikp,ephmat)
7 use modmain
8 use modphonon
9 implicit none
10 ! arguments
11 integer, intent(in) :: iqp,ikp
12 complex(8), intent(out) :: ephmat(nstsv,nstsv,nbph)
13 ! local variables
14 integer isym,lspl,iq,ik,iv(3)
15 integer recl,n,nstsv_,nbph_
16 real(8) vql_(3),vkl_(3),t1
17 if (iqp <= nqpt) then
18 ! q-point is in the reduced set
19  iq=iqp
20  ik=ikp
21 else
22 ! q-point is not in the reduced set
23  call findqpt(vql(:,iqp),isym,iq)
24  lspl=lsplsymc(isym)
25  call i3mtv(symlat(:,:,lspl),ivk(:,ikp),iv)
26  iv(:)=modulo(iv(:),ngridk(:))
27  ik=ivkiknr(iv(1),iv(2),iv(3))
28 end if
29 ! find the record length
30 inquire(iolength=recl) vql_,vkl_,nstsv_,nbph_,ephmat
31 ! record number
32 n=(iq-1)*nkptnr+ik
33 !$OMP CRITICAL(u240)
34 open(240,file='EPHMAT.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
35 read(240,rec=n) vql_,vkl_,nstsv_,nbph_,ephmat
36 close(240)
37 !$OMP END CRITICAL(u240)
38 t1=abs(vql(1,iq)-vql_(1))+abs(vql(2,iq)-vql_(2))+abs(vql(3,iq)-vql_(3))
39 if (t1 > epslat) then
40  write(*,*)
41  write(*,'("Error(getephmat): differing vectors for q-point ",I8)') iq
42  write(*,'(" current : ",3G18.10)') vql(:,iq)
43  write(*,'(" EPHMAT.OUT : ",3G18.10)') vql_
44  write(*,*)
45  stop
46 end if
47 t1=abs(vkl(1,ik)-vkl_(1))+abs(vkl(2,ik)-vkl_(2))+abs(vkl(3,ik)-vkl_(3))
48 if (t1 > epslat) then
49  write(*,*)
50  write(*,'("Error(getephmat): differing vectors for k-point ",I8)') ik
51  write(*,'(" current : ",3G18.10)') vkl(:,ik)
52  write(*,'(" EPHMAT.OUT : ",3G18.10)') vkl_
53  write(*,*)
54  stop
55 end if
56 if (nstsv /= nstsv_) then
57  write(*,*)
58  write(*,'("Error(getephmat): differing nstsv for q-point, k-point ",2I8)') &
59  iq,ik
60  write(*,'(" current : ",I8)') nstsv
61  write(*,'(" EPHMAT.OUT : ",I8)') nstsv_
62  write(*,*)
63  stop
64 end if
65 if (nbph /= nbph_) then
66  write(*,*)
67  write(*,'("Error(getephmat): differing nbph for q-point, k-point ",2I8)') &
68  iq,ik
69  write(*,'(" current : ",I8)') nbph
70  write(*,'(" EPHMAT.OUT : ",I8)') nbph_
71  write(*,*)
72  stop
73 end if
74 end subroutine
75 
integer nqpt
Definition: modmain.f90:525
integer, dimension(3, 3, 48) symlat
Definition: modmain.f90:344
integer nkptnr
Definition: modmain.f90:463
pure subroutine i3mtv(a, x, y)
Definition: i3mtv.f90:10
integer, dimension(maxsymcrys) lsplsymc
Definition: modmain.f90:364
subroutine findqpt(vpl, isym, iq)
Definition: findqpt.f90:7
integer, dimension(3) ngridk
Definition: modmain.f90:448
real(8), dimension(:,:), allocatable vql
Definition: modmain.f90:545
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
real(8) epslat
Definition: modmain.f90:24
subroutine getephmat(iqp, ikp, ephmat)
Definition: getephmat.f90:7
integer, dimension(:,:), allocatable ivk
Definition: modmain.f90:465
integer, dimension(:,:,:), allocatable ivkiknr
Definition: modmain.f90:469