File: readdyn.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (55 lines) | stat: -rw-r--r-- 1,325 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55

! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine readdyn(dynq)
use modmain
use modphonon
implicit none
! arguments
complex(8), intent(out) :: dynq(nbph,nbph,nqpt)
! local variables
logical exist
integer iq,is,js,ia,ja
integer ip,jp,i,j
real(8) a,b
character(256) fext
do iq=1,nqpt
  i=0
  do is=1,nspecies
    do ia=1,natoms(is)
      do ip=1,3
        i=i+1
        call phfext(iq,is,ia,ip,fext)
        inquire(file='DYN'//trim(fext),exist=exist)
        if (.not.exist) then
          write(*,*)
          write(*,'("Error(readdyn): file not found :")')
          write(*,'(A)') ' DYN'//trim(fext)
          write(*,*)
          stop
        end if
        open(50,file='DYN'//trim(fext),status='OLD',form='FORMATTED')
        j=0
        do js=1,nspecies
          do ja=1,natoms(js)
            do jp=1,3
              j=j+1
              read(50,*) a,b
              dynq(i,j,iq)=cmplx(a,b,8)
            end do
          end do
        end do
        close(50)
      end do
! end loops over atoms and species
    end do
  end do
! symmetrise the dynamical matrix
  call dynsym(vql(:,iq),dynq(:,:,iq))
! end loop over q-vectors
end do
return
end subroutine