The Elk Code
 
Loading...
Searching...
No Matches
hermite.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: hermite
8! !INTERFACE:
9real(8) function hermite(n,x)
10! !INPUT/OUTPUT PARAMETERS:
11! n : order of Hermite polynomial (in,integer)
12! x : real argument (in,real)
13! !DESCRIPTION:
14! Returns the $n$th Hermite polynomial. The recurrence relation
15! $$ H_i(x)=2xH_{i-1}(x)-2nH_{i-2}(x), $$
16! with $H_0=1$ and $H_1=2x$, is used. This procedure is numerically stable
17! and accurate to near machine precision for $n\le 20$.
18!
19! !REVISION HISTORY:
20! Created April 2003 (JKD)
21!EOP
22!BOC
23implicit none
24! arguments
25integer, intent(in) :: n
26real(8), intent(in) :: x
27! local variables
28integer i
29real(8) h1,h2,ht
30! fast return if possible
31if (n == 0) then
32 hermite=1.d0
33 return
34else if (n == 1) then
35 hermite=2.d0*x
36 return
37else if (n == 2) then
38 hermite=4.d0*x**2-2.d0
39 return
40end if
41if (n < 0) then
42 write(*,*)
43 write(*,'("Error(hermite): n < 0 : ",I8)') n
44 write(*,*)
45 stop
46end if
47if (n > 20) then
48 write(*,*)
49 write(*,'("Error(hermite): n out of range : ",I8)') n
50 write(*,*)
51 stop
52end if
53if (abs(x) > 1.d15) then
54 write(*,*)
55 write(*,'("Error(hermite): x out of range : ",G18.10)') x
56 write(*,*)
57 stop
58end if
59h1=2.d0*x
60h2=1.d0
61do i=2,n
62 ht=2.d0*(x*h1-dble(i-1)*h2)
63 h2=h1
64 h1=ht
65end do
66hermite=h1
67end function
68!EOC
69
real(8) function hermite(n, x)
Definition hermite.f90:10