The Elk Code
 
Loading...
Searching...
No Matches
stheta.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: stheta
8! !INTERFACE:
9real(8) function stheta(stype,x)
10! !INPUT/OUTPUT PARAMETERS:
11! stype : smearing type (in,integer)
12! x : real argument (in,real)
13! !DESCRIPTION:
14! Returns the Heaviside step function corresponding to the smooth
15! approximation to the Dirac delta function:
16! $$ \tilde\Theta(x)=\int_{-\infty}^x dt\,\tilde\delta(t). $$
17! See function {\tt sdelta} for details.
18!
19! !REVISION HISTORY:
20! Created April 2003 (JKD)
21!EOP
22!BOC
23implicit none
24! arguments
25integer, intent(in) :: stype
26real(8), intent(in) :: x
27! external functions
28real(8), external :: stheta_mp,stheta_fd,stheta_sq,stheta_lr
29select case(stype)
30case(0)
31 stheta=stheta_mp(0,x)
32case(1)
33 stheta=stheta_mp(1,x)
34case(2)
35 stheta=stheta_mp(2,x)
36case(3)
38case(4)
40case(5)
42case default
43 write(*,*)
44 write(*,'("Error(stheta): stype not defined : ",I8)') stype
45 write(*,*)
46 stop
47end select
48end function
49!EOC
50
real(8) function stheta(stype, x)
Definition stheta.f90:10
elemental real(8) function stheta_fd(x)
Definition stheta_fd.f90:10
elemental real(8) function stheta_lr(x)
Definition stheta_lr.f90:7
real(8) function stheta_mp(n, x)
Definition stheta_mp.f90:10
elemental real(8) function stheta_sq(x)
Definition stheta_sq.f90:10