The Elk Code
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 
6 subroutine phlwidth
7 use modmain
8 use modphonon
9 implicit none
10 ! local variables
11 integer iq,i,j,iv
12 real(8) gmin,gmax,t1
13 ! allocatable arrays
14 real(8), allocatable :: wq(:),gmr(:,:,:),gq(:,:),gp(:,:)
15 complex(8), allocatable :: dq(:,:),ev(:,:),b(:,:)
16 complex(8), allocatable :: gmq(:,:,:),gmp(:,:)
17 ! initialise universal variables
18 call init0
19 call init2
20 call initph
21 allocate(wq(nbph),gmr(nbph,nbph,nqptnr),gq(nbph,nqpt),gp(nbph,npp1d))
22 allocate(dq(nbph,nbph),ev(nbph,nbph),b(nbph,nbph))
23 allocate(gmq(nbph,nbph,nqpt),gmp(nbph,nbph))
24 ! read in the phonon linewidths for each q-point
25 call readgamma(gq)
26 ! loop over phonon q-points
27 do 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)
39 end do
40 ! Fourier transform the gamma matrices to real-space
41 call dynqtor(gmq,gmr)
42 ! generate a set of q-point vectors along a path in the Brillouin zone
44 gmin=1.d8
45 gmax=0.d0
46 ! compute the linewidths along the path
47 do 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
63 end do
64 t1=(gmax-gmin)*0.5d0
65 gmin=gmin-t1
66 gmax=gmax+t1
67 ! output the vertex location lines
68 open(50,file='PHLWLINES.OUT',form='FORMATTED')
69 do iv=1,nvp1d
70  write(50,'(2G18.10)') dvp1d(iv),gmin
71  write(50,'(2G18.10)') dvp1d(iv),gmax
72  write(50,*)
73 end do
74 close(50)
75 ! output the phonon linewidth dispersion
76 open(50,file='PHLWIDTH.OUT',form='FORMATTED')
77 do i=1,nbph
78  do iq=ip01d,npp1d
79  write(50,'(2G18.10)') dpp1d(iq),gp(i,iq)
80  end do
81  write(50,*)
82 end do
83 close(50)
84 write(*,*)
85 write(*,'("Info(phlwidth):")')
86 write(*,'(" phonon linewidth dispersion written to PHLWIDTH.OUT")')
87 write(*,'(" vertex location lines written to PHLWLINES.OUT")')
88 deallocate(wq,gmr,gq,gp,dq,ev,b,gmq,gmp)
89 end subroutine
90 
real(8), dimension(:,:,:), allocatable dynr
Definition: modphonon.f90:29
integer nqpt
Definition: modmain.f90:525
real(8), dimension(:), allocatable dpp1d
Definition: modmain.f90:1126
subroutine readgamma(gq)
Definition: readgamma.f90:7
complex(8), parameter zone
Definition: modmain.f90:1240
subroutine dynrtoq(vpl, dr, dq)
Definition: dynrtoq.f90:7
real(8), dimension(:), allocatable dvp1d
Definition: modmain.f90:1122
subroutine dynev(dq, wq, ev)
Definition: dynev.f90:7
complex(8), dimension(:,:,:), allocatable dynq
Definition: modphonon.f90:27
integer nbph
Definition: modphonon.f90:13
subroutine init2
Definition: init2.f90:7
subroutine initph
Definition: initph.f90:7
subroutine plotpt1d(cvec, nv, np, vvl, vpl, dv, dp)
Definition: plotpt1d.f90:10
integer ip01d
Definition: modmain.f90:1118
real(8), dimension(3, 3) bvec
Definition: modmain.f90:16
integer nqptnr
Definition: modmain.f90:527
complex(8), parameter zzero
Definition: modmain.f90:1240
real(8), dimension(:,:), allocatable vvlp1d
Definition: modmain.f90:1120
subroutine dynqtor(dq, dr)
Definition: dynqtor.f90:7
subroutine init0
Definition: init0.f90:10
integer nvp1d
Definition: modmain.f90:1114
integer npp1d
Definition: modmain.f90:1116
subroutine phlwidth
Definition: phlwidth.f90:7
real(8), dimension(:,:), allocatable vplp1d
Definition: modmain.f90:1124
subroutine dynevs(ev, dq, wq)
Definition: dynevs.f90:7