File: mcmillan.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (65 lines) | stat: -rw-r--r-- 1,532 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65

! Copyright (C) 2011 J. K. Dewhurst, A. Sanna, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine mcmillan(w,a2f,lambda,wlog,wrms,tc)
use modmain
use modphonon
implicit none
! arguments
real(8), intent(in) :: w(nwplot),a2f(nwplot)
real(8), intent(out) :: lambda,wlog,wrms,tc
! local variables
integer iw
real(8) l1,l2,f1,f2,t1
! allocatable arrays
real(8), allocatable :: f(:)
! external functions
real(8) fintgt
external fintgt
allocate(f(nwplot))
! compute the total lambda
do iw=1,nwplot
  if (w(iw).gt.1.d-8) then
    f(iw)=a2f(iw)/w(iw)
  else
    f(iw)=0.d0
  end if
end do
t1=fintgt(-3,nwplot,w,f)
lambda=2.d0*t1
! compute the logarithmic average frequency
do iw=1,nwplot
  if (w(iw).gt.1.d-8) then
    f(iw)=a2f(iw)*log(w(iw))/w(iw)
  else
    f(iw)=0.d0
  end if
end do
t1=fintgt(-3,nwplot,w,f)
t1=(2.d0/lambda)*t1
wlog=exp(t1)
! compute < w^2 >^(1/2)
do iw=1,nwplot
  if (w(iw).gt.1.d-8) then
    f(iw)=a2f(iw)*w(iw)
  else
    f(iw)=0.d0
  end if
end do
t1=fintgt(-3,nwplot,w,f)
t1=(2.d0/lambda)*t1
wrms=sqrt(abs(t1))
! compute McMillan-Allen-Dynes superconducting critical temperature
t1=(-1.04d0*(1.d0+lambda))/(lambda-mustar-0.62d0*lambda*mustar)
tc=(wlog/(1.2d0*kboltz))*exp(t1)
l1=2.46d0*(1.d0+3.8d0*mustar)
l2=1.82d0*(1.d0+6.3d0*mustar)*(wrms/wlog)
f1=(1.d0+(lambda/l1)**(3.d0/2.d0))**(1.d0/3.d0)
f2=1.d0+(wrms/wlog-1.d0)*(lambda**2)/(lambda**2+l2**2)
tc=tc*f1*f2
deallocate(f)
return
end subroutine