File: sdelta.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 (113 lines) | stat: -rw-r--r-- 2,680 bytes parent folder | download
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
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113

! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU Lesser General Public
! License. See the file COPYING for license details.

!BOP
! !ROUTINE: sdelta
! !INTERFACE:
real(8) function sdelta(stype,x)
! !INPUT/OUTPUT PARAMETERS:
!   stype : smearing type (in,integer)
!   x     : real argument (in,real)
! !DESCRIPTION:
!   Returns a normalised smooth approximation to the Dirac delta function. These
!   functions are defined such that
!   $$ \int\tilde{\delta}(x)dx=1. $$
!   The effective width, $w$, of the delta function may be varied by using the
!   normalising transformation
!   $$ \tilde{\delta}_w(x)\equiv\frac{\tilde{\delta}(x/w)}{w}. $$
!   Currently implimented are:
!   \begin{list}{}{\itemsep -2pt}
!    \item[0.] Gaussian
!    \item[1.] Methfessel-Paxton order 1
!    \item[2.] Methfessel-Paxton order 2
!    \item[3.] Fermi-Dirac
!    \item[4.] Square-wave impulse
!    \item[5.] Lorentzian
!   \end{list}
!   See routines {\tt stheta}, {\tt sdelta\_mp}, {\tt sdelta\_fd} and
!   {\tt sdelta\_sq}.
!
! !REVISION HISTORY:
!   Created April 2003 (JKD)
!EOP
!BOC
implicit none
! arguments
integer, intent(in) :: stype
real(8), intent(in) :: x
! external functions
real(8) sdelta_mp,sdelta_fd,sdelta_sq,sdelta_lr
external sdelta_mp,sdelta_fd,sdelta_sq,sdelta_lr
select case(stype)
case(0)
  sdelta=sdelta_mp(0,x)
  return
case(1)
  sdelta=sdelta_mp(1,x)
  return
case(2)
  sdelta=sdelta_mp(2,x)
  return
case(3)
  sdelta=sdelta_fd(x)
  return
case(4)
  sdelta=sdelta_sq(x)
case(5)
  sdelta=sdelta_lr(x)
case default
  write(*,*)
  write(*,'("Error(sdelta): sytpe not defined : ",I8)') stype
  write(*,*)
  stop
end select
end function
!EOC

!BOP
! !ROUTINE: getsdata
! !INTERFACE:
subroutine getsdata(stype,sdescr)
! !INPUT/OUTPUT PARAMETERS:
!   stype  : smearing type (in,integer)
!   sdescr : smearing scheme description (out,character(*))
! !DESCRIPTION:
!   Returns a description of the smearing scheme as string {\tt sdescr} up to
!   256 characters long.
!
! !REVISION HISTORY:
!   Created April 2003 (JKD)
!EOP
!BOC
implicit none
! arguments
integer, intent(in) :: stype
character(*), intent(out) :: sdescr
select case(stype)
case(0)
  sdescr='Gaussian'
  return
case(1)
  sdescr='Methfessel-Paxton order 1, Phys. Rev. B 40, 3616 (1989)'
  return
case(2)
  sdescr='Methfessel-Paxton order 2, Phys. Rev. B 40, 3616 (1989)'
  return
case(3)
  sdescr='Fermi-Dirac'
  return
case(4)
  sdescr='Square-wave impulse'
case(5)
  sdescr='Lorentzian'
case default
  write(*,*)
  write(*,'("Error(getsdata): sytpe not defined : ",I8)') stype
  write(*,*)
  stop
end select
end subroutine
!EOC