File: phdisp.f90

package info (click to toggle)
elkcode 10.8.1-1
  • links: PTS, VCS
  • area: main
  • in suites: experimental
  • size: 10,672 kB
  • sloc: f90: 52,747; perl: 964; makefile: 352; sh: 296; python: 105; ansic: 67
file content (60 lines) | stat: -rw-r--r-- 1,655 bytes parent folder | download | duplicates (2)
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
56
57
58
59
60

! 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 phdisp
use modmain
use modphonon
implicit none
! local variables
integer iq,i,iv
real(8) wmin,wmax,t1
! allocatable arrays
real(8), allocatable :: wq(:,:)
complex(8), allocatable :: dq(:,:),ev(:,:)
! initialise universal variables
call init0
call init2
call initph
! allocate local arrays
allocate(wq(nbph,npp1d),dq(nbph,nbph),ev(nbph,nbph))
! generate a set of q-point vectors along a path in the Brillouin zone
call plotpt1d(bvec,nvp1d,npp1d,vvlp1d,vplp1d,dvp1d,dpp1d)
! compute the phonon frequencies along the path
do iq=ip01d,npp1d
! compute the dynamical matrix at this particular q-point
  call dynrtoq(vplp1d(:,iq),dynr,dq)
! find the phonon frequencies and eigenvectors
  call dynev(dq,wq(:,iq),ev)
end do
! find the minimum and maximum phonon frequencies
wmin=minval(wq(1,:))
wmax=maxval(wq(nbph,:))
t1=(wmax-wmin)*0.5d0
wmin=wmin-t1
wmax=wmax+t1
! output the vertex location lines
open(50,file='PHDLINES.OUT',form='FORMATTED')
do iv=1,nvp1d
  write(50,'(2G18.10)') dvp1d(iv),wmin
  write(50,'(2G18.10)') dvp1d(iv),wmax
  write(50,*)
end do
close(50)
! output the phonon dispersion
open(50,file='PHDISP.OUT',form='FORMATTED')
do i=1,nbph
  do iq=ip01d,npp1d
    write(50,'(2G18.10)') dpp1d(iq),wq(i,iq)
  end do
  write(50,*)
end do
close(50)
write(*,*)
write(*,'("Info(phdisp):")')
write(*,'(" phonon dispersion written to PHDISP.OUT")')
write(*,'(" vertex location lines written to PHDLINES.OUT")')
deallocate(dq,ev)
end subroutine