File: hermite.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (67 lines) | stat: -rw-r--r-- 1,390 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67

! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU Lesser General Public
! License. See the file COPYING for license details.

!BOP
! !ROUTINE: hermite
! !INTERFACE:
real(8) function hermite(n,x)
! !INPUT/OUTPUT PARAMETERS:
!   n : order of Hermite polynomial (in,integer)
!   x : real argument (in,real)
! !DESCRIPTION:
!   Returns the $n$th Hermite polynomial. The recurrence relation
!   $$ H_i(x)=2xH_{i-1}(x)-2nH_{i-2}(x), $$
!   with $H_0=1$ and $H_1=2x$, is used. This procedure is numerically stable
!   and accurate to near machine precision for $n\le 20$.
!
! !REVISION HISTORY:
!   Created April 2003 (JKD)
!EOP
!BOC
implicit none
! arguments
integer, intent(in) :: n
real(8), intent(in) :: x
! local variables
integer i
real(8) h1,h2,ht
if (n.lt.0) then
  write(*,*)
  write(*,'("Error(hermite): n < 0 : ",I8)') n
  write(*,*)
  stop
end if
if (n.gt.20) then
  write(*,*)
  write(*,'("Error(hermite): n out of range : ",I8)') n
  write(*,*)
  stop
end if
if (abs(x).gt.1.d15) then
  write(*,*)
  write(*,'("Error(hermite): x out of range : ",G18.10)') x
  write(*,*)
  stop
end if
if (n.eq.0) then
  hermite=1.d0
  return
end if
if (n.eq.1) then
  hermite=2.d0*x
  return
end if
h1=2.d0*x
h2=1.d0
do i=2,n
  ht=2.d0*x*h1-2.d0*dble(i-1)*h2
  h2=h1
  h1=ht
end do
hermite=h1
return
end function
!EOC