The Elk Code
 
Loading...
Searching...
No Matches
phdos.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 phdos
7use modmain
8use modphonon
9use modtest
10implicit none
11! local variables
12integer iq,i,iw
13integer i1,i2,i3
14real(8) wmin,wmax,wd,dw
15real(8) tmax,temp(ntemp),s(ntemp)
16real(8) v(3),t1,t2
17! allocatable arrays
18real(8), allocatable :: wq(:),w(:),gw(:),f(:)
19complex(8), allocatable :: dq(:,:),ev(:,:)
20! external functions
21real(8), external :: splint
22! initialise universal variables
23call init0
24call init2
25call initph
26allocate(wq(nbph),w(nwplot),gw(nwplot),f(nwplot))
27allocate(dq(nbph,nbph),ev(nbph,nbph))
28! find the minimum and maximum frequencies
29wmin=0.d0
30wmax=0.d0
31do iq=1,nqpt
32 call dynev(dynq(:,:,iq),wq,ev)
33 wmin=min(wmin,wq(1))
34 wmax=max(wmax,wq(nbph))
35end do
36t1=(wmax-wmin)*0.1d0
37wmin=wmin-t1
38wmax=wmax+t1
39wd=wmax-wmin
40if (wd < 1.d-8) wd=1.d0
41dw=wd/dble(nwplot)
42! generate energy grid
43do iw=1,nwplot
44 w(iw)=dw*dble(iw-1)+wmin
45end do
46gw(:)=0.d0
47do i1=0,ngrkf-1
48 v(1)=dble(i1)/dble(ngrkf)
49 do i2=0,ngrkf-1
50 v(2)=dble(i2)/dble(ngrkf)
51 do i3=0,ngrkf-1
52 v(3)=dble(i3)/dble(ngrkf)
53! compute the dynamical matrix at this particular q-point
54 call dynrtoq(v,dynr,dq)
55! find the phonon frequencies
56 call dynev(dq,wq,ev)
57 do i=1,nbph
58 t1=(wq(i)-wmin)/dw+1.d0
59 iw=nint(t1)
60 if ((iw >= 1).and.(iw <= nwplot)) then
61 gw(iw)=gw(iw)+1.d0
62 end if
63 end do
64 end do
65 end do
66end do
67t1=1.d0/(dw*dble(ngrkf)**3)
68gw(:)=t1*gw(:)
69! smooth phonon DOS if required
70if (nswplot > 0) call fsmooth(nswplot,nwplot,gw)
71! write phonon DOS to file
72open(50,file='PHDOS.OUT',form='FORMATTED')
73do iw=1,nwplot
74 write(50,'(2G18.10)') w(iw),gw(iw)
75end do
76close(50)
77write(*,*)
78write(*,'("Info(phdos):")')
79write(*,'(" phonon density of states written to PHDOS.OUT")')
80!-------------------------------------------!
81! thermodynamic properties from DOS !
82!-------------------------------------------!
83! maximum temperature
84tmax=wmax/kboltz
85! temperature grid
86do i=1,ntemp
87 temp(i)=tmax*dble(i)/dble(ntemp)
88end do
89open(50,file='THERMO.OUT',form='FORMATTED')
90write(50,*)
91write(50,'("(All quantities are per unit cell)")')
92! zero point energy
93do iw=1,nwplot
94 f(iw)=gw(iw)*w(iw)
95end do
96t1=splint(nwplot,w,f)
97t1=0.5d0*t1
98write(50,*)
99write(50,'("Zero-point energy : ",G18.10)') t1
100! vibrational energy
101write(50,*)
102write(50,'("Vibrational energy vs. temperature :")')
103do i=1,ntemp
104 do iw=1,nwplot
105 t1=w(iw)/(2.d0*kboltz*temp(i))
106 if (t1 > 0.d0) then
107 f(iw)=gw(iw)*w(iw)*cosh(t1)/sinh(t1)
108 else
109 f(iw)=0.d0
110 end if
111 end do
112 t1=splint(nwplot,w,f)
113 t1=0.5d0*t1
114 write(50,'(2G18.10)') temp(i),t1
115 s(i)=t1
116end do
117! free energy
118write(50,*)
119write(50,'("Free energy vs. temperature :")')
120do i=1,ntemp
121 do iw=1,nwplot
122 t1=2.d0*sinh(w(iw)/(2.d0*kboltz*temp(i)))
123 if (t1 > 0.d0) then
124 f(iw)=gw(iw)*log(t1)
125 else
126 f(iw)=0.d0
127 end if
128 end do
129 t1=splint(nwplot,w,f)
130 t1=kboltz*temp(i)*t1
131 write(50,'(2G18.10)') temp(i),t1
132! compute entropy from S = (F-E)/T
133 s(i)=(s(i)-t1)/temp(i)
134end do
135! entropy
136write(50,*)
137write(50,'("Entropy vs. temperature :")')
138do i=1,ntemp
139 write(50,'(2G18.10)') temp(i),s(i)
140end do
141! heat capacity
142write(50,*)
143write(50,'("Heat capacity vs. temperature :")')
144do i=1,ntemp
145 do iw=1,nwplot
146 t1=w(iw)/(kboltz*temp(i))
147 t2=exp(t1)-1.d0
148 if (abs(t2) > 1.d-14) then
149 f(iw)=gw(iw)*(t1**2)*(t2+1.d0)/t2**2
150 else
151 f(iw)=0.d0
152 end if
153 end do
154 t1=splint(nwplot,w,f)
155 t1=kboltz*t1
156 write(50,'(2G18.10)') temp(i),t1
157end do
158close(50)
159write(*,'(" thermodynamic properties written to THERMO.OUT")')
160! write phonon DOS to test file
161call writetest(210,'phonon DOS',nv=nwplot,tol=1.d-2,rva=gw)
162deallocate(wq,w,gw,f,dq,ev)
163end subroutine
164
subroutine dynev(dq, wq, ev)
Definition dynev.f90:7
subroutine dynrtoq(vpl, dr, dq)
Definition dynrtoq.f90:7
pure subroutine fsmooth(m, n, f)
Definition fsmooth.f90:10
subroutine init0
Definition init0.f90:10
subroutine init2
Definition init2.f90:7
subroutine initph
Definition initph.f90:7
integer ngrkf
Definition modmain.f90:1072
integer nwplot
Definition modmain.f90:1070
integer nswplot
Definition modmain.f90:1074
integer nqpt
Definition modmain.f90:525
real(8), parameter kboltz
Definition modmain.f90:1262
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 writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
subroutine phdos
Definition phdos.f90:7
real(8) function splint(n, x, f)
Definition splint.f90:7