The Elk Code
Loading...
Searching...
No Matches
polynm.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: polynm
8
! !INTERFACE:
9
pure
real(8)
function
polynm
(m,np,xa,ya,x)
10
! !INPUT/OUTPUT PARAMETERS:
11
! m : order of derivative (in,integer)
12
! np : number of points to fit (in,integer)
13
! ip : point at which m'th derivative is to be evaluated (in,integer)
14
! xa : abscissa array (in,real(np))
15
! ya : ordinate array (in,real(np))
16
! x : evaluation abscissa (in,real)
17
! !DESCRIPTION:
18
! Fits a polynomial of order $n_p-1$ to a set of $n_p$ points. If $m\ge 0$ the
19
! function returns the $m$th derviative of the polynomial at $x$, while for
20
! $m<0$ the integral of the polynomial from the first point in the array to
21
! $x$ is returned.
22
!
23
! !REVISION HISTORY:
24
! Created October 2002 (JKD)
25
! Simplified, January 2025 (JKD)
26
!EOP
27
!BOC
28
implicit none
29
! argmuments
30
integer
,
intent(in)
:: m,np
31
real
(8),
intent(in)
:: xa(np),ya(np),x
32
! local variables
33
integer
i,j,k
34
real
(8) x0,t1
35
! automatic arrays
36
real
(8) c(np)
37
if
((np < 1).or.(m >= np))
then
38
polynm
=0.d0
39
return
40
end if
41
! find the polynomial coefficients in divided differences form
42
c(1:np)=ya(1:np)
43
do
i=2,np
44
do
j=np,i,-1
45
c(j)=(c(j)-c(j-1))/(xa(j)-xa(j+1-i))
46
end do
47
end do
48
! special case m = 0
49
if
(m == 0)
then
50
polynm
=c(1)
51
t1=1.d0
52
do
i=2,np
53
t1=t1*(x-xa(i-1))
54
polynm
=
polynm
+c(i)*t1
55
end do
56
return
57
end if
58
x0=xa(1)
59
! convert to standard form
60
do
j=1,np-1
61
do
i=1,np-j
62
k=np-i
63
c(k)=c(k)+(x0-xa(k-j+1))*c(k+1)
64
end do
65
end do
66
if
(m > 0)
then
67
! take the m'th derivative
68
do
j=1,m
69
do
i=m+1,np
70
c(i)=c(i)*dble(i-j)
71
end do
72
end do
73
polynm
=c(np)
74
t1=x-x0
75
do
i=np-1,m+1,-1
76
polynm
=
polynm
*t1+c(i)
77
end do
78
else
79
! find the integral
80
polynm
=c(np)/dble(np)
81
t1=x-x0
82
do
i=np-1,1,-1
83
polynm
=
polynm
*t1+c(i)/dble(i)
84
end do
85
polynm
=
polynm
*t1
86
end if
87
end function
88
!EOC
89
polynm
pure real(8) function polynm(m, np, xa, ya, x)
Definition
polynm.f90:10
polynm.f90
Generated by
1.9.8