The Elk Code
readgamma.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2008 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 readgamma(gq)
7 use modmain
8 use modphonon
9 implicit none
10 ! arguments
11 real(8), intent(out) :: gq(nbph,nqpt)
12 ! local variables
13 integer iq,i
14 integer natmtot_,nqpt_,iq_,i_
15 real(8) vql_(3),vqc_(3),t1
16 open(50,file='GAMMAQ.OUT',form='FORMATTED',status='OLD')
17 read(50,*)
18 read(50,*) natmtot_
19 if (natmtot /= natmtot_) then
20  write(*,*)
21  write(*,'("Error(readgamma): differing natmtot")')
22  write(*,'(" current : ",I4)') natmtot
23  write(*,'(" GAMMAQ.OUT : ",I4)') natmtot_
24  write(*,*)
25  stop
26 end if
27 read(50,*) nqpt_
28 if (nqpt /= nqpt_) then
29  write(*,*)
30  write(*,'("Error(readgamma): differing nqpt")')
31  write(*,'(" current : ",I6)') nqpt
32  write(*,'(" GAMMAQ.OUT : ",I6)') nqpt_
33  write(*,*)
34  stop
35 end if
36 read(50,*)
37 do iq=1,nqpt
38  read(50,*) iq_
39  if (iq /= iq_) then
40  write(*,*)
41  write(*,'("Error(readgamma): incorrect q-point index in GAMMAQ.OUT for &
42  &q-point ",I6)') iq
43  write(*,*)
44  stop
45  end if
46  read(50,*) vql_
47  t1=sum(abs(vql(:,iq)-vql_(:)))
48  if (t1 > epslat) then
49  write(*,*)
50  write(*,'("Error(readgamma): differing q-vectors in lattice coordinates &
51  &for q-point ",I6)') iq
52  write(*,'(" current : ",3G18.10)') vql(:,iq)
53  write(*,'(" GAMMAQ.OUT : ",3G18.10)') vql_
54  write(*,*)
55  stop
56  end if
57  read(50,*) vqc_
58  t1=sum(abs(vqc(:,iq)-vqc_(:)))
59  if (t1 > epslat) then
60  write(*,*)
61  write(*,'("Error(readgamma): differing q-vectors in Cartesian coordinates &
62  &for q-point ",I6)') iq
63  write(*,'(" current : ",3G18.10)') vqc(:,iq)
64  write(*,'(" GAMMAQ.OUT : ",3G18.10)') vqc_
65  write(*,*)
66  stop
67  end if
68  do i=1,nbph
69  read(50,*) i_,gq(i,iq)
70  if (i /= i_) then
71  write(*,*)
72  write(*,'("Error(readgamma): incorrect mode index in GAMMAQ.OUT for &
73  &q-point ",I6)') iq
74  write(*,*)
75  stop
76  end if
77  end do
78  read(50,*)
79 end do
80 close(50)
81 end subroutine
82 
subroutine readgamma(gq)
Definition: readgamma.f90:7
real(8), dimension(:,:), allocatable vqc
Definition: modmain.f90:547
real(8), dimension(:,:), allocatable vql
Definition: modmain.f90:545
real(8) epslat
Definition: modmain.f90:24
integer natmtot
Definition: modmain.f90:40