The Elk Code
 
Loading...
Searching...
No Matches
sdelta.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 Lesser General Public
4! License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: sdelta
8! !INTERFACE:
9real(8) function sdelta(stype,x)
10! !INPUT/OUTPUT PARAMETERS:
11! stype : smearing type (in,integer)
12! x : real argument (in,real)
13! !DESCRIPTION:
14! Returns a normalised smooth approximation to the Dirac delta function. These
15! functions are defined such that
16! $$ \int\tilde{\delta}(x)dx=1. $$
17! The effective width, $w$, of the delta function may be varied by using the
18! normalising transformation
19! $$ \tilde{\delta}_w(x)\equiv\frac{\tilde{\delta}(x/w)}{w}. $$
20! Currently implimented are:
21! \begin{list}{}{\itemsep -2pt}
22! \item[0.] Gaussian
23! \item[1.] Methfessel-Paxton order 1
24! \item[2.] Methfessel-Paxton order 2
25! \item[3.] Fermi-Dirac
26! \item[4.] Square-wave impulse
27! \item[5.] Lorentzian
28! \end{list}
29! See routines {\tt stheta}, {\tt sdelta\_mp}, {\tt sdelta\_fd} and
30! {\tt sdelta\_sq}.
31!
32! !REVISION HISTORY:
33! Created April 2003 (JKD)
34!EOP
35!BOC
36implicit none
37! arguments
38integer, intent(in) :: stype
39real(8), intent(in) :: x
40! external functions
41real(8), external :: sdelta_mp,sdelta_fd,sdelta_sq,sdelta_lr
42select case(stype)
43case(0)
44 sdelta=sdelta_mp(0,x)
45case(1)
46 sdelta=sdelta_mp(1,x)
47case(2)
48 sdelta=sdelta_mp(2,x)
49case(3)
51case(4)
53case(5)
55case default
56 write(*,*)
57 write(*,'("Error(sdelta): stype not defined : ",I8)') stype
58 write(*,*)
59 stop
60end select
61end function
62!EOC
63
64!BOP
65! !ROUTINE: getsdata
66! !INTERFACE:
67subroutine getsdata(stype,sdescr)
68! !INPUT/OUTPUT PARAMETERS:
69! stype : smearing type (in,integer)
70! sdescr : smearing scheme description (out,character(*))
71! !DESCRIPTION:
72! Returns a description of the smearing scheme as string {\tt sdescr} up to
73! 256 characters long.
74!
75! !REVISION HISTORY:
76! Created April 2003 (JKD)
77!EOP
78!BOC
79implicit none
80! arguments
81integer, intent(in) :: stype
82character(*), intent(out) :: sdescr
83select case(stype)
84case(0)
85 sdescr='Gaussian'
86case(1)
87 sdescr='Methfessel-Paxton order 1, Phys. Rev. B 40, 3616 (1989)'
88case(2)
89 sdescr='Methfessel-Paxton order 2, Phys. Rev. B 40, 3616 (1989)'
90case(3)
91 sdescr='Fermi-Dirac'
92case(4)
93 sdescr='Square-wave impulse'
94case(5)
95 sdescr='Lorentzian'
96case default
97 write(*,*)
98 write(*,'("Error(getsdata): stype not defined : ",I8)') stype
99 write(*,*)
100 stop
101end select
102end subroutine
103!EOC
104
real(8) function sdelta(stype, x)
Definition sdelta.f90:10
subroutine getsdata(stype, sdescr)
Definition sdelta.f90:68
elemental real(8) function sdelta_fd(x)
Definition sdelta_fd.f90:10
elemental real(8) function sdelta_lr(x)
Definition sdelta_lr.f90:7
real(8) function sdelta_mp(n, x)
Definition sdelta_mp.f90:10
elemental real(8) function sdelta_sq(x)
Definition sdelta_sq.f90:10