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