-
Notifications
You must be signed in to change notification settings - Fork 0
/
spliq.f90
239 lines (175 loc) · 6.52 KB
/
spliq.f90
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
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
!> spliq integrates a cubic spline (generated by
!> splift, smoo, etc.) on the intervals (xlo,xup(i)), where xup
!> is a sequence of upper limits on the intervals of integration.
!> the only restrictions on xlo and xup(*) are
!> xlo < xup(1),
!> xup(i) .le. xup(i+1) for each i .
!> endpoints beyond the span of abscissas are allowed.
!> the spline over the interval (x(i),x(i+1)) is regarded
!> as a cubic polynomial expanded about x(i) and is integrated analytically.
!>
!> The code is based on, and compatible with, the 1978 subroutine
!> written by M. K. Gordon at Sandia.
!>
!> \author Jose Luis Martins
!> \version 5.805
!> \date 9 May 2018
!> \copyright GNU Public License v2
subroutine spliq(x,y,yp,ypp,n,xlo,xup,nup,ans,ierr)
! Bug nup = 1 squashed, 1 October 2020.
! copyright José Luís Martins, INESC-MN.
implicit none
integer, parameter :: REAL64 = selected_real_kind(12)
integer, parameter :: iowrite = 6
! input
integer, intent(in) :: n !< number of data points that define the spline. n > 1
real(REAL64), intent(in) :: x(n) !< array of abscissas (in increasing order) that define the spline.
real(REAL64), intent(in) :: y(n) !< array of ordinates that define the spline.
real(REAL64), intent(in) :: yp(n) !< array of first derivatives of the spline (at the x(i)).
real(REAL64), intent(in) :: ypp(n) !< array of second derivatives that define the spline.
real(REAL64), intent(in) :: xlo !< left endpoint of integration intervals
integer, intent(in) :: nup !< number of right endpoints
real(REAL64), intent(in) :: xup(nup) !< array of right endpoints in ascending order
! output
real(REAL64), intent(out) :: ans(nup) !< array of integral values, that is, ans(i) = integral from xlo to xup(i)
integer, intent(out) :: ierr !< status code. 1: all normal;
!< 2: improper input - n < 4 or nup < 1;
!< 3: improper input - abscissas not in strictly ascending order
!< 4: improper input - right endpoints xup not in ascending order
!< 5: improper input - xlo > xup(1)
!< 6: integration successful but with extrapolation
! local variables
integer :: i ! i is current index into x array.
real(REAL64) :: hlo, hlo2, hi, hi2, hi3
real(REAL64) :: hup, hup2, hup3, hup4
real(REAL64) :: hsum, hdiff
real(REAL64) :: sumx ! accumulates integral values
real(REAL64) :: sum0, sum1, sum2, sum3
real(REAL64) :: psum0, psum1, psum2, psum3
integer :: kmin
logical :: lexit
! counter
integer :: j, k, m
! constants
real(REAL64), parameter :: ZERO = 0.0_REAL64, UM = 1.0_REAL64
! ck for improper input
ierr = 0
if(n < 4 .or. nup < 1) then
ierr = 2
endif
if(ierr == 0) then
do k = 1,n-1
if(x(k) >= x(k+1)) then
ierr = 3
exit
endif
enddo
endif
if(ierr == 0) then
if(nup /= 1) then
do k = 2,nup
if(xup(k-1) > xup(k)) then
ierr = 4
exit
endif
enddo
endif
endif
if(ierr == 0) then
if(xlo > xup(1)) then
ierr = 5
endif
endif
if(ierr == 0) then
if(xlo < x(1) .or. xup(nup) > x(n)) then
ierr = 6
else
ierr = 1
endif
! ocate xlo in interval (x(i),x(i+1))
i = n-1
do j = 1,n-2
if(xlo < x(j+1)) then
i = j
exit
endif
enddo
hlo = xlo-x(i)
hlo2 = hlo*hlo
hi = x(i+1)-x(i)
hi2 = hi*hi
j = 1
lexit = .FALSE.
do m = 1,nup
if(xup(m) > x(i+1) .and. xlo < x(n-1)) then
j = m
lexit = .TRUE.
exit
else
! compute special cases of xup in interval with xlo
hup = xup(m)-x(i)
hsum = hup+hlo
hdiff = hup-hlo
hup2 = hup*hup
sumx = (ypp(i+1)-ypp(i))*hsum*hdiff*(hup2+hlo2)/(24*hi)
sumx = sumx + ypp(i)*hdiff*(hup2+hlo*hup+hlo2)/6
sumx = sumx + yp(i)*hdiff*hsum/2
sumx = sumx + y(i)*hdiff
ans(m) = sumx
j = m
endif
enddo
! xup(nup) outside interval with xlo
if(lexit) then
! compute integral between xlo and x(i+1) as four terms in taylor
! polynomial and advance i to i+1
hdiff = hi-hlo
hsum = hi+hlo
sum0 = y(i)*hdiff
sum1 = yp(i)*hdiff*hsum
sum2 = ypp(i)*hdiff*(hi2+hi*hlo+hlo2)
sum3 = (ypp(i+1)-ypp(i))*hdiff*hsum*(hi2+hlo2)/hi
i = i+1
! locate each xup(m) in interval (x(i),x(i+1))
do m = j,nup
kmin = i
do k = kmin,n
if(xup(m) < x(k+1) .or. k == n-1) then
i = k
exit
endif
! augment integral between abscissas to include interval
! (x(k),x(k+1))
hi = x(k+1)-x(k)
hi2 = hi*hi
hi3 = hi2*hi
sum0 = sum0 + y(k)*hi
sum1 = sum1 + yp(k)*hi2
sum2 = sum2 + ypp(k)*hi3
sum3 = sum3 + (ypp(k+1)-ypp(k))*hi3
enddo
! integral between x(i) and xup(m) is zero
if(xup(m) /= x(i)) then
! compute integral between x(i) and xup(m) and evaluate
! taylor polynomial in reverse order
hup = xup(m)-x(i)
hup2 = hup*hup
hup3 = hup2*hup
hup4 = hup3*hup
hi = x(i+1)-x(i)
psum0 = y(i)*hup
psum1 = yp(i)*hup2
psum2 = ypp(i)*hup3
psum3 = (ypp(i+1)-ypp(i))*hup4/hi
sumx = (sum3+psum3)/24 + (sum2+psum2)/6
sumx = sumx + (sum1+psum1)/2
sumx = sumx + (sum0+psum0)
else
sumx = ((sum3/24 + sum2/6) + sum1/2) + sum0
endif
ans(m) = sumx
enddo
endif
endif
return
end subroutine spliq