The Elk Code
Loading...
Searching...
No Matches
phdisp.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2002-2005 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
phdisp
7
use
modmain
8
use
modphonon
9
implicit none
10
! local variables
11
integer
iq,i,iv
12
real
(8) wmin,wmax,t1
13
! allocatable arrays
14
real
(8),
allocatable
:: wq(:,:)
15
complex(8)
,
allocatable
:: dq(:,:),ev(:,:)
16
! initialise universal variables
17
call
init0
18
call
init2
19
call
initph
20
! allocate local arrays
21
allocate
(wq(
nbph
,
npp1d
),dq(
nbph
,
nbph
),ev(
nbph
,
nbph
))
22
! generate a set of q-point vectors along a path in the Brillouin zone
23
call
plotpt1d
(
bvec
,
nvp1d
,
npp1d
,
vvlp1d
,
vplp1d
,
dvp1d
,
dpp1d
)
24
! compute the phonon frequencies along the path
25
do
iq=
ip01d
,
npp1d
26
! compute the dynamical matrix at this particular q-point
27
call
dynrtoq
(
vplp1d
(:,iq),
dynr
,dq)
28
! find the phonon frequencies and eigenvectors
29
call
dynev
(dq,wq(:,iq),ev)
30
end do
31
! find the minimum and maximum phonon frequencies
32
wmin=minval(wq(1,:))
33
wmax=maxval(wq(
nbph
,:))
34
t1=(wmax-wmin)*0.5d0
35
wmin=wmin-t1
36
wmax=wmax+t1
37
! output the vertex location lines
38
open
(50,file=
'PHDLINES.OUT'
,form=
'FORMATTED'
)
39
do
iv=1,
nvp1d
40
write
(50,
'(2G18.10)'
)
dvp1d
(iv),wmin
41
write
(50,
'(2G18.10)'
)
dvp1d
(iv),wmax
42
write
(50,*)
43
end do
44
close
(50)
45
! output the phonon dispersion
46
open
(50,file=
'PHDISP.OUT'
,form=
'FORMATTED'
)
47
do
i=1,
nbph
48
do
iq=
ip01d
,
npp1d
49
write
(50,
'(2G18.10)'
)
dpp1d
(iq),wq(i,iq)
50
end do
51
write
(50,*)
52
end do
53
close
(50)
54
write
(*,*)
55
write
(*,
'("Info(phdisp):")'
)
56
write
(*,
'(" phonon dispersion written to PHDISP.OUT")'
)
57
write
(*,
'(" vertex location lines written to PHDLINES.OUT")'
)
58
deallocate
(dq,ev)
59
end subroutine
60
dynev
subroutine dynev(dq, wq, ev)
Definition
dynev.f90:7
dynrtoq
subroutine dynrtoq(vpl, dr, dq)
Definition
dynrtoq.f90:7
init0
subroutine init0
Definition
init0.f90:10
init2
subroutine init2
Definition
init2.f90:7
initph
subroutine initph
Definition
initph.f90:7
modmain
Definition
modmain.f90:6
modmain::npp1d
integer npp1d
Definition
modmain.f90:1113
modmain::nvp1d
integer nvp1d
Definition
modmain.f90:1111
modmain::bvec
real(8), dimension(3, 3) bvec
Definition
modmain.f90:16
modmain::vvlp1d
real(8), dimension(:,:), allocatable vvlp1d
Definition
modmain.f90:1117
modmain::ip01d
integer ip01d
Definition
modmain.f90:1115
modmain::vplp1d
real(8), dimension(:,:), allocatable vplp1d
Definition
modmain.f90:1121
modmain::dvp1d
real(8), dimension(:), allocatable dvp1d
Definition
modmain.f90:1119
modmain::dpp1d
real(8), dimension(:), allocatable dpp1d
Definition
modmain.f90:1123
modphonon
Definition
modphonon.f90:6
modphonon::nbph
integer nbph
Definition
modphonon.f90:13
modphonon::dynr
real(8), dimension(:,:,:), allocatable dynr
Definition
modphonon.f90:29
phdisp
subroutine phdisp
Definition
phdisp.f90:7
plotpt1d
subroutine plotpt1d(cvec, nv, np, vvl, vpl, dv, dp)
Definition
plotpt1d.f90:10
phdisp.f90
Generated by
1.9.8