The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine mcmillan(w,a2f,lambda,wlog,wrms,tc)
7use modmain
8use modphonon
9implicit none
10! arguments
11real(8), intent(in) :: w(nwplot),a2f(nwplot)
12real(8), intent(out) :: lambda,wlog,wrms,tc
13! local variables
14integer iw
15real(8) l1,l2,f1,f2,t1
16! allocatable arrays
17real(8), allocatable :: f(:)
18! external functions
19real(8), external :: splint
20allocate(f(nwplot))
21! compute the total lambda
22do 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
28end do
29t1=splint(nwplot,w,f)
30lambda=2.d0*t1
31! compute the logarithmic average frequency
32do 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
38end do
39t1=splint(nwplot,w,f)
40t1=(2.d0/lambda)*t1
41wlog=exp(t1)
42! compute ⟨w²⟩¹⸍²
43do 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
49end do
50t1=splint(nwplot,w,f)
51t1=(2.d0/lambda)*t1
52wrms=sqrt(abs(t1))
53! compute McMillan-Allen-Dynes superconducting critical temperature
54t1=(-1.04d0*(1.d0+lambda))/(lambda-mustar-0.62d0*lambda*mustar)
55tc=(wlog/(1.2d0*kboltz))*exp(t1)
56l1=2.46d0*(1.d0+3.8d0*mustar)
57l2=1.82d0*(1.d0+6.3d0*mustar)*(wrms/wlog)
58f1=(1.d0+(lambda/l1)**(3.d0/2.d0))**(1.d0/3.d0)
59f2=1.d0+(wrms/wlog-1.d0)*(lambda**2)/(lambda**2+l2**2)
60tc=tc*f1*f2
61deallocate(f)
62end subroutine
63
subroutine mcmillan(w, a2f, lambda, wlog, wrms, tc)
Definition mcmillan.f90:7
real(8), parameter kboltz
Definition modmain.f90:1262
real(8) mustar
Definition modphonon.f90:23
real(8) function splint(n, x, f)
Definition splint.f90:7