The Elk Code
 
Loading...
Searching...
No Matches
phlwidth.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
6subroutine phlwidth
7use modmain
8use modphonon
9implicit none
10! local variables
11integer iq,i,j,iv
12real(8) gmin,gmax,t1
13! allocatable arrays
14real(8), allocatable :: wq(:),gmr(:,:,:),gq(:,:),gp(:,:)
15complex(8), allocatable :: dq(:,:),ev(:,:),b(:,:)
16complex(8), allocatable :: gmq(:,:,:),gmp(:,:)
17! initialise universal variables
18call init0
19call init2
20call initph
21allocate(wq(nbph),gmr(nbph,nbph,nqptnr),gq(nbph,nqpt),gp(nbph,npp1d))
22allocate(dq(nbph,nbph),ev(nbph,nbph),b(nbph,nbph))
23allocate(gmq(nbph,nbph,nqpt),gmp(nbph,nbph))
24! read in the phonon linewidths for each q-point
25call readgamma(gq)
26! loop over phonon q-points
27do iq=1,nqpt
28! find the eigenvalues and eigenvectors of the dynamical matrix
29 call dynev(dynq(:,:,iq),wq,ev)
30! construct a complex matrix from the phonon eigenvectors such that its
31! eigenvalues squared are the phonon linewidths
32 do i=1,nbph
33 t1=sqrt(abs(gq(i,iq)))
34 do j=1,nbph
35 b(i,j)=t1*conjg(ev(j,i))
36 end do
37 end do
38 call zgemm('N','N',nbph,nbph,nbph,zone,ev,nbph,b,nbph,zzero,gmq(:,:,iq),nbph)
39end do
40! Fourier transform the gamma matrices to real-space
41call dynqtor(gmq,gmr)
42! generate a set of q-point vectors along a path in the Brillouin zone
44gmin=1.d8
45gmax=0.d0
46! compute the linewidths along the path
47do iq=ip01d,npp1d
48! compute the dynamical matrix at this particular q-point
49 call dynrtoq(vplp1d(:,iq),dynr,dq)
50! find the phonon eigenvalues and eigenvectors
51 call dynev(dq,wq,ev)
52! compute the gamma matrix at this particular q-point
53 call dynrtoq(vplp1d(:,iq),gmr,gmp)
54! diagonalise the gamma matrix simultaneously with the dynamical matrix
55 call dynevs(ev,gmp,gp(:,iq))
56! square the eigenvalues to recover the linewidths
57 gp(:,iq)=gp(:,iq)**2
58 do i=1,nbph
59 t1=gp(i,iq)
60 if (t1 < gmin) gmin=t1
61 if (t1 > gmax) gmax=t1
62 end do
63end do
64t1=(gmax-gmin)*0.5d0
65gmin=gmin-t1
66gmax=gmax+t1
67! output the vertex location lines
68open(50,file='PHLWLINES.OUT',form='FORMATTED')
69do iv=1,nvp1d
70 write(50,'(2G18.10)') dvp1d(iv),gmin
71 write(50,'(2G18.10)') dvp1d(iv),gmax
72 write(50,*)
73end do
74close(50)
75! output the phonon linewidth dispersion
76open(50,file='PHLWIDTH.OUT',form='FORMATTED')
77do i=1,nbph
78 do iq=ip01d,npp1d
79 write(50,'(2G18.10)') dpp1d(iq),gp(i,iq)
80 end do
81 write(50,*)
82end do
83close(50)
84write(*,*)
85write(*,'("Info(phlwidth):")')
86write(*,'(" phonon linewidth dispersion written to PHLWIDTH.OUT")')
87write(*,'(" vertex location lines written to PHLWLINES.OUT")')
88deallocate(wq,gmr,gq,gp,dq,ev,b,gmq,gmp)
89end subroutine
90
subroutine dynev(dq, wq, ev)
Definition dynev.f90:7
subroutine dynevs(ev, dq, wq)
Definition dynevs.f90:7
subroutine dynqtor(dq, dr)
Definition dynqtor.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
complex(8), parameter zzero
Definition modmain.f90:1238
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 nqptnr
Definition modmain.f90:527
complex(8), parameter zone
Definition modmain.f90:1238
integer ip01d
Definition modmain.f90:1115
integer nqpt
Definition modmain.f90:525
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
complex(8), dimension(:,:,:), allocatable dynq
Definition modphonon.f90:27
subroutine phlwidth
Definition phlwidth.f90:7
subroutine plotpt1d(cvec, nv, np, vvl, vpl, dv, dp)
Definition plotpt1d.f90:10
subroutine readgamma(gq)
Definition readgamma.f90:7