-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathgrid_interp.f90
327 lines (237 loc) · 7.41 KB
/
grid_interp.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
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
!> Interpolates from one grid to another grid using Lagrange interpolation.
!> Returns also the derivatives of the function on the new grid
!> and an estimate of the respective errors.
!> Based on D. B. Hunter, The Computer Journal, 3, 270 (1961)
!> with the small difference algorithm of Numerical Recipes Eq. 3.1.5
!>
!> \author Jose Luis Martins
!> \version 6.0.6 of atom
!> \date 15 July 2021.
!> \copyright GNU Public License v3
subroutine grid_interp(n, nd, nin, xin, fin, nout, xout, fout, dymax, ifail, mxdout)
implicit none
integer, parameter :: REAL64 = selected_real_kind(12)
integer, parameter :: iowrite = 6 ! default tape for writing
! input
integer, intent(in) :: mxdout !< dimension of points in the output (new) grid
integer, intent(in) :: n !< order of Lagrange interpolation
integer, intent(in) :: nd !< order of derivatives
integer, intent(in) :: nin !< number of points in the input (old) grid
real(REAL64), intent(in) :: xin(nin) !< points on the input grid
real(REAL64), intent(in) :: fin(nin) !< f(x(i))
integer, intent(in) :: nout !< number of points in the output (new) grid
real(REAL64), intent(in) :: xout(mxdout) !< points on the output (new) grid
! output
real(REAL64), intent(out) :: fout(mxdout,0:nd) !< fout(:,0): value of the function on the new grid; fout(:,n): n-th derivative
real(REAL64), intent(out) :: dymax(0:nd) !< estimate of the accuracy (from n-1 interpolation)
integer, intent(out) :: ifail !< ifail =/= 0 indicates failure of the process.
! allocatable arrays
real(REAL64), allocatable :: xs(:) ! local old shifted grid
real(REAL64), allocatable :: ys(:) ! function on old shifted grid
real(REAL64), allocatable :: y(:) ! f(0): interpolated value at x=0; f(n): interpolated derivative at x=0.
real(REAL64), allocatable :: dy(:) ! error estimate (last correction)
integer, allocatable :: isin(:), isout(:) ! indexation of grids
! local variables
logical :: linup, lindw
logical :: loutup, loutdw
real(REAL64) :: xdif
integer :: near
integer :: nn
real(REAL64) :: xgmin, xgmax
integer :: ninf, nsup
! constants
real(REAL64), parameter :: ZERO = 0.0_REAL64
real(REAL64), parameter :: EPS = 1.0E-9_REAL64
! counters
integer :: i, j
! check grid sizes and if original grid is ordered
ifail = 0
if(nin < 2 .or. nin < n+1) then
write(iowrite,*) ' grid_interp, nin = ', nin, ' n = ',n
ifail = 1
return
endif
if(n < 1) then
write(iowrite,*) ' grid_interp, n = ', n
ifail = 2
return
endif
if(nout > mxdout) then
write(iowrite,*) ' grid_interp, nout = ', nout, ' mxdout = ', mxdout
ifail = 3
return
endif
if(n > 10) then
write(iowrite,*) ' grid_interp, n = ', n
write(iowrite,*) ' Do you know what you are doing? '
endif
xgmin = xin(1)
xgmax = xin(1)
do i = 2,nin
if(xgmin > xin(i)) xgmin = xin(i)
if(xgmax < xin(i)) xgmax = xin(i)
enddo
xdif = EPS*(xgmax - xgmin)
do j = 1,nout
if(xout(j) < xgmin-xdif .or. xout(j) > xgmax+xdif) then
write(iowrite,*) ' grid_interp, extrapolation not allowed'
write(iowrite,*) ' xout(j),xgmin,xgmax = ', xout(j),xgmin,xgmax
ifail = 4
return
endif
enddo
! check if the input grid is ordered. Close points are not allowed
linup = .TRUE.
do i = 1,nin-1
if(xin(i+1) < xin(i)+xdif) then
linup = .FALSE.
exit
endif
enddo
lindw = .TRUE.
do i = 1,nin-1
if(xin(i+1)+xdif > xin(i)) then
lindw = .FALSE.
exit
endif
enddo
! check if the output grid is ordered
loutup = .TRUE.
do i = 1,nout-1
if(xout(i+1) < xout(i)) then
loutup = .FALSE.
exit
endif
enddo
loutdw = .TRUE.
do i = 1,nout-1
if(xout(i+1) > xout(i)) then
loutdw = .FALSE.
exit
endif
enddo
! allocate arrays
allocate(xs(0:n),ys(0:n))
allocate(y(0:nd),dy(0:nd))
do i = 0,nd
dymax(i) = ZERO
enddo
if( (linup .or. lindw) .and. (loutup .or. loutdw) ) then
if( (linup .and. loutup) .or. (lindw .and. loutdw) ) then
near = 1
else
near = nin
endif
do j = 1,nout
! find nearest point below
if(linup .and. loutup) then
nn = near
do i = near,nin
if(xin(i) > xout(j)) exit
nn = i
enddo
near = nn
elseif(lindw .and. loutup) then
nn = near
do i = near,1,-1
nn = i
if(xin(i) > xout(j)) exit
enddo
near = nn
elseif(linup .and. loutdw) then
nn = near
do i = near,1,-1
nn = i
if(xin(i) < xout(j)) exit
enddo
near = nn
elseif(lindw .and. loutdw) then
nn = near
do i = near,nin
if(xin(i) < xout(j)) exit
nn = i
enddo
near = nn
endif
ninf = near - n/2
if(ninf < 1) ninf = 1
nsup = ninf + n
if(nsup > nin) then
nsup = nin
ninf = nsup - n
endif
if(linup) then
do i = 0,n
xs(i) = xin(ninf+i) - xout(j)
ys(i) = fin(ninf+i)
enddo
else
do i = 0,n
xs(i) = xin(nsup-i) - xout(j)
ys(i) = fin(nsup-i)
enddo
endif
call poly_interp(y, dy, xs, ys, n, nd)
do i = 0,nd
fout(j,i) = y(i)
if(dymax(i) < dy(i)) dymax(i) = dy(i)
enddo
enddo
else
! one or both of the grids is not ordered
allocate(isin(nin),isout(nout))
if(linup) then
do j = 1,nin
isin(j) = j
enddo
else
call sort(nin,xin,isin)
endif
if(loutup) then
do j = 1,nout
isout(j) = j
enddo
else
call sort(nout,xout,isout)
endif
! check if the input grid has close points
do i = 1,nin-1
if( xin(isin(i+1)) < xin(isin(i))+xdif ) then
write(iowrite,*) ' grid_interp, duplicate poins in input grid'
write(6,*) 'xin(',isin(i+1),') = ',xin(isin(i+1))
write(6,*) 'xin(',isin(i),') = ',xin(isin(i))
ifail = 5
return
endif
enddo
near = 1
do j = 1,nout
! find nearest point below
nn = near
do i = near,nin
if(xin(isin(i)) > xout(isout(j))) exit
nn = i
enddo
near = nn
ninf = near - n/2
if(ninf < 1) ninf = 1
nsup = ninf + n
if(nsup > nin) then
nsup = nin
ninf = nsup - n
endif
do i = 0,n
xs(i) = xin(isin(ninf+i)) - xout(isout(j))
ys(i) = fin(isin(ninf+i))
enddo
call poly_interp(y, dy, xs, ys, n, nd)
do i = 0,nd
fout(isout(j),i) = y(i)
if(dymax(i) < dy(i)) dymax(i) = dy(i)
enddo
enddo
endif
deallocate(xs,ys)
deallocate(y,dy)
return
end subroutine grid_interp