The Elk Code
 
Loading...
Searching...
No Matches
sdelta_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: sdelta_mp
8! !INTERFACE:
9real(8) function sdelta_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 Dirac delta function of order $N$
15! given by Methfessel and Paxton, {\it Phys. Rev. B} {\bf 40}, 3616 (1989),
16! $$ \tilde\delta(x)=\sum_{i=0}^N \frac{(-1)^i}{i!4^n\sqrt\pi} H_{2i}(x)
17! e^{-x^2},$$
18! where $H_j$ is the $j$th-order Hermite polynomial. This function has the
19! property
20! $$ \int_{-\infty}^{\infty}\tilde\delta(x)P(x)=P(0), $$
21! where $P(x)$ is any polynomial of degree $2N+1$ or less. The case $N=0$
22! corresponds to Gaussian smearing. This procedure is numerically stable
23! and accurate to near machine precision for $N\le 10$.
24!
25! !REVISION HISTORY:
26! Created April 2003 (JKD)
27!EOP
28!BOC
29implicit none
30! arguments
31integer, intent(in) :: n
32real(8), intent(in) :: x
33! local variables
34integer i
35real(8), parameter :: sqpi=1.7724538509055160273d0
36real(8) sm,t0,t1
37! external functions
38real(8), external :: factn,hermite
39if (n == 0) then
40 sdelta_mp=exp(-x**2)/sqpi
41 return
42end if
43if (n < 0) then
44 write(*,*)
45 write(*,'("Error(sdelta_mp): n < 0 : ",I8)') n
46 write(*,*)
47 stop
48end if
49if (n > 10) then
50 write(*,*)
51 write(*,'("Error(sdelta_mp): n out of range : ",I8)') n
52 write(*,*)
53 stop
54end if
55if (abs(x) > 12.d0) then
56 sdelta_mp=0.d0
57 return
58end if
59t0=exp(-x**2)/sqpi
60sm=t0
61do i=1,n
62 t1=t0/(factn(i)*dble(4**i))
63 if (mod(i,2) /= 0) t1=-t1
64 sm=sm+t1*hermite(2*i,x)
65end do
66sdelta_mp=sm
67end function
68!EOC
69
elemental real(8) function factn(n)
Definition factn.f90:7
real(8) function hermite(n, x)
Definition hermite.f90:10
real(8) function sdelta_mp(n, x)
Definition sdelta_mp.f90:10