The Elk Code
stheta_mp.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_mp
8 ! !INTERFACE:
9 real(8) function stheta_mp(n,x)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! n : order (in,integer)
12 ! x : real argument (in,real)
13 ! !DESCRIPTION:
14 ! Returns the smooth approximation to the Heaviside step function of order
15 ! $N$ given by Methfessel and Paxton, {\it Phys. Rev. B} {\bf 40}, 3616
16 ! (1989),
17 ! $$ \tilde\Theta(x)=1-S_N(x) $$
18 ! where
19 ! \begin{align*}
20 ! S_N(x)&=S_0(x)+\sum_{i=1}^N \frac{(-1)^i}{i!4^n\sqrt\pi} H_{2i-1}(x)
21 ! e^{-x^2},\\
22 ! S_0(x)&=\frac{1}{2}(1-{\rm erf}(x))
23 ! \end{align*}
24 ! and $H_j$ is the $j$th-order Hermite polynomial. This procedure is
25 ! numerically stable and accurate to near machine precision for $N\le 10$.
26 !
27 ! !REVISION HISTORY:
28 ! Created April 2003 (JKD)
29 !EOP
30 !BOC
31 implicit none
32 ! arguments
33 integer, intent(in) :: n
34 real(8), intent(in) :: x
35 ! local variables
36 integer i
37 real(8), parameter :: sqpi=1.7724538509055160273d0
38 real(8) sm,t0,t1
39 ! external functions
40 real(8), external :: factn,hermite,erf
41 if (n == 0) then
42  stheta_mp=0.5d0*(1.d0+erf(x))
43  return
44 end if
45 if (n < 0) then
46  write(*,*)
47  write(*,'("Error(stheta_mp): n < 0 : ",I8)') n
48  write(*,*)
49  stop
50 end if
51 if (n > 10) then
52  write(*,*)
53  write(*,'("Error(stheta_mp): n out of range : ",I8)') n
54  write(*,*)
55  stop
56 end if
57 if (x < -12.d0) then
58  stheta_mp=0.d0
59  return
60 end if
61 if (x > 12.d0) then
62  stheta_mp=1.d0
63  return
64 end if
65 t0=-exp(-x**2)/sqpi
66 sm=0.5d0*(1.d0+erf(x))
67 do i=1,n
68  t1=t0/(factn(i)*dble(4**i))
69  if (mod(i,2) /= 0) t1=-t1
70  sm=sm+t1*hermite(2*i-1,x)
71 end do
72 stheta_mp=sm
73 end function
74 !EOC
75 
real(8) function stheta_mp(n, x)
Definition: stheta_mp.f90:10
real(8) function hermite(n, x)
Definition: hermite.f90:10
elemental real(8) function erf(x)
Definition: erf.f90:9
elemental real(8) function factn(n)
Definition: factn.f90:7