The Elk Code
mcmillan.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2011 J. K. Dewhurst, A. Sanna, S. Sharma and E. K. U. Gross.
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 mcmillan(w,a2f,lambda,wlog,wrms,tc)
7 use modmain
8 use modphonon
9 implicit none
10 ! arguments
11 real(8), intent(in) :: w(nwplot),a2f(nwplot)
12 real(8), intent(out) :: lambda,wlog,wrms,tc
13 ! local variables
14 integer iw
15 real(8) l1,l2,f1,f2,t1
16 ! allocatable arrays
17 real(8), allocatable :: f(:)
18 ! external functions
19 real(8), external :: splint
20 allocate(f(nwplot))
21 ! compute the total lambda
22 do iw=1,nwplot
23  if (w(iw) > 1.d-8) then
24  f(iw)=a2f(iw)/w(iw)
25  else
26  f(iw)=0.d0
27  end if
28 end do
29 t1=splint(nwplot,w,f)
30 lambda=2.d0*t1
31 ! compute the logarithmic average frequency
32 do iw=1,nwplot
33  if (w(iw) > 1.d-8) then
34  f(iw)=a2f(iw)*log(w(iw))/w(iw)
35  else
36  f(iw)=0.d0
37  end if
38 end do
39 t1=splint(nwplot,w,f)
40 t1=(2.d0/lambda)*t1
41 wlog=exp(t1)
42 ! compute ⟨w²⟩¹⸍²
43 do iw=1,nwplot
44  if (w(iw) > 1.d-8) then
45  f(iw)=a2f(iw)*w(iw)
46  else
47  f(iw)=0.d0
48  end if
49 end do
50 t1=splint(nwplot,w,f)
51 t1=(2.d0/lambda)*t1
52 wrms=sqrt(abs(t1))
53 ! compute McMillan-Allen-Dynes superconducting critical temperature
54 t1=(-1.04d0*(1.d0+lambda))/(lambda-mustar-0.62d0*lambda*mustar)
55 tc=(wlog/(1.2d0*kboltz))*exp(t1)
56 l1=2.46d0*(1.d0+3.8d0*mustar)
57 l2=1.82d0*(1.d0+6.3d0*mustar)*(wrms/wlog)
58 f1=(1.d0+(lambda/l1)**(3.d0/2.d0))**(1.d0/3.d0)
59 f2=1.d0+(wrms/wlog-1.d0)*(lambda**2)/(lambda**2+l2**2)
60 tc=tc*f1*f2
61 deallocate(f)
62 end subroutine
63 
subroutine mcmillan(w, a2f, lambda, wlog, wrms, tc)
Definition: mcmillan.f90:7
real(8), parameter kboltz
Definition: modmain.f90:1263
real(8) mustar
Definition: modphonon.f90:23