1 |
!! |
2 |
!! Copyright (c) 2006 The University of Notre Dame. All Rights Reserved. |
3 |
!! |
4 |
!! The University of Notre Dame grants you ("Licensee") a |
5 |
!! non-exclusive, royalty free, license to use, modify and |
6 |
!! redistribute this software in source and binary code form, provided |
7 |
!! that the following conditions are met: |
8 |
!! |
9 |
!! 1. Redistributions of source code must retain the above copyright |
10 |
!! notice, this list of conditions and the following disclaimer. |
11 |
!! |
12 |
!! 2. Redistributions in binary form must reproduce the above copyright |
13 |
!! notice, this list of conditions and the following disclaimer in the |
14 |
!! documentation and/or other materials provided with the |
15 |
!! distribution. |
16 |
!! |
17 |
!! This software is provided "AS IS," without a warranty of any |
18 |
!! kind. All express or implied conditions, representations and |
19 |
!! warranties, including any implied warranty of merchantability, |
20 |
!! fitness for a particular purpose or non-infringement, are hereby |
21 |
!! excluded. The University of Notre Dame and its licensors shall not |
22 |
!! be liable for any damages suffered by licensee as a result of |
23 |
!! using, modifying or distributing the software or its |
24 |
!! derivatives. In no event will the University of Notre Dame or its |
25 |
!! licensors be liable for any lost revenue, profit or data, or for |
26 |
!! direct, indirect, special, consequential, incidental or punitive |
27 |
!! damages, however caused and regardless of the theory of liability, |
28 |
!! arising out of the use of or inability to use software, even if the |
29 |
!! University of Notre Dame has been advised of the possibility of |
30 |
!! such damages. |
31 |
!! |
32 |
!! SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your |
33 |
!! research, please cite the appropriate papers when you publish your |
34 |
!! work. Good starting points are: |
35 |
!! |
36 |
!! [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). |
37 |
!! [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). |
38 |
!! [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008). |
39 |
!! [4] Vardeman & Gezelter, in progress (2009). |
40 |
!! |
41 |
!! |
42 |
!! interpolation.F90 |
43 |
!! |
44 |
!! Created by Charles F. Vardeman II on 03 Apr 2006. |
45 |
!! |
46 |
!! PURPOSE: Generic Spline interpolation routines. |
47 |
!! |
48 |
!! @author Charles F. Vardeman II |
49 |
!! @version $Id$ |
50 |
|
51 |
|
52 |
module interpolation |
53 |
use definitions |
54 |
use status |
55 |
implicit none |
56 |
PRIVATE |
57 |
|
58 |
type, public :: cubicSpline |
59 |
logical :: isUniform = .false. |
60 |
integer :: n = 0 |
61 |
real(kind=dp) :: dx_i |
62 |
real (kind=dp), pointer,dimension(:) :: x => null() |
63 |
real (kind=dp), pointer,dimension(:) :: y => null() |
64 |
real (kind=dp), pointer,dimension(:) :: b => null() |
65 |
real (kind=dp), pointer,dimension(:) :: c => null() |
66 |
real (kind=dp), pointer,dimension(:) :: d => null() |
67 |
end type cubicSpline |
68 |
|
69 |
public :: newSpline |
70 |
public :: deleteSpline |
71 |
public :: lookupSpline |
72 |
public :: lookupUniformSpline |
73 |
public :: lookupNonuniformSpline |
74 |
public :: lookupUniformSpline1d |
75 |
|
76 |
contains |
77 |
|
78 |
|
79 |
subroutine newSpline(cs, x, y, isUniform) |
80 |
|
81 |
implicit none |
82 |
|
83 |
type (cubicSpline), intent(inout) :: cs |
84 |
real( kind = DP ), intent(in) :: x(:), y(:) |
85 |
real( kind = DP ) :: fp1, fpn, p |
86 |
REAL( KIND = DP), DIMENSION(size(x)-1) :: diff_y, H |
87 |
|
88 |
logical, intent(in) :: isUniform |
89 |
integer :: i, alloc_error, n, k |
90 |
|
91 |
alloc_error = 0 |
92 |
|
93 |
if (cs%n .ne. 0) then |
94 |
call handleWarning("interpolation::newSpline", & |
95 |
"cubicSpline struct was already created") |
96 |
call deleteSpline(cs) |
97 |
end if |
98 |
|
99 |
! make sure the sizes match |
100 |
|
101 |
n = size(x) |
102 |
|
103 |
if ( size(y) .ne. size(x) ) then |
104 |
call handleError("interpolation::newSpline", & |
105 |
"Array size mismatch") |
106 |
end if |
107 |
|
108 |
cs%n = n |
109 |
cs%isUniform = isUniform |
110 |
|
111 |
allocate(cs%x(n), stat=alloc_error) |
112 |
if(alloc_error .ne. 0) then |
113 |
call handleError("interpolation::newSpline", & |
114 |
"Error in allocating storage for x") |
115 |
endif |
116 |
|
117 |
allocate(cs%y(n), stat=alloc_error) |
118 |
if(alloc_error .ne. 0) then |
119 |
call handleError("interpolation::newSpline", & |
120 |
"Error in allocating storage for y") |
121 |
endif |
122 |
|
123 |
allocate(cs%b(n), stat=alloc_error) |
124 |
if(alloc_error .ne. 0) then |
125 |
call handleError("interpolation::newSpline", & |
126 |
"Error in allocating storage for b") |
127 |
endif |
128 |
|
129 |
allocate(cs%c(n), stat=alloc_error) |
130 |
if(alloc_error .ne. 0) then |
131 |
call handleError("interpolation::newSpline", & |
132 |
"Error in allocating storage for c") |
133 |
endif |
134 |
|
135 |
allocate(cs%d(n), stat=alloc_error) |
136 |
if(alloc_error .ne. 0) then |
137 |
call handleError("interpolation::newSpline", & |
138 |
"Error in allocating storage for d") |
139 |
endif |
140 |
|
141 |
! make sure we are monotinically increasing in x: |
142 |
|
143 |
h = diff(x) |
144 |
if (any(h <= 0)) then |
145 |
call handleError("interpolation::newSpline", & |
146 |
"Negative dx interval found") |
147 |
end if |
148 |
|
149 |
! load x and y values into the cubicSpline structure: |
150 |
|
151 |
do i = 1, n |
152 |
cs%x(i) = x(i) |
153 |
cs%y(i) = y(i) |
154 |
end do |
155 |
|
156 |
! Calculate coefficients for the tridiagonal system: store |
157 |
! sub-diagonal in B, diagonal in D, difference quotient in C. |
158 |
|
159 |
cs%b(1:n-1) = h |
160 |
diff_y = diff(y) |
161 |
cs%c(1:n-1) = diff_y / h |
162 |
|
163 |
if (n == 2) then |
164 |
! Assume the derivatives at both endpoints are zero |
165 |
! another assumption could be made to have a linear interpolant |
166 |
! between the two points. In that case, the b coefficients |
167 |
! below would be diff_y(1)/h(1) and the c and d coefficients would |
168 |
! both be zero. |
169 |
cs%b(1) = 0.0_dp |
170 |
cs%c(1) = -3.0_dp * (diff_y(1)/h(1))**2 |
171 |
cs%d(1) = -2.0_dp * (diff_y(1)/h(1))**3 |
172 |
cs%b(2) = cs%b(1) |
173 |
cs%c(2) = 0.0_dp |
174 |
cs%d(2) = 0.0_dp |
175 |
cs%dx_i = 1.0_dp / h(1) |
176 |
return |
177 |
end if |
178 |
|
179 |
cs%d(1) = 2.0_dp * cs%b(1) |
180 |
do i = 2, n-1 |
181 |
cs%d(i) = 2.0_dp * (cs%b(i) + cs%b(i-1)) |
182 |
end do |
183 |
cs%d(n) = 2.0_dp * cs%b(n-1) |
184 |
|
185 |
! Calculate estimates for the end slopes using polynomials |
186 |
! that interpolate the data nearest the end. |
187 |
|
188 |
fp1 = cs%c(1) - cs%b(1)*(cs%c(2) - cs%c(1))/(cs%b(1) + cs%b(2)) |
189 |
if (n > 3) then |
190 |
fp1 = fp1 + cs%b(1)*((cs%b(1) + cs%b(2))*(cs%c(3) - cs%c(2))/ & |
191 |
(cs%b(2) + cs%b(3)) - cs%c(2) + cs%c(1))/(x(4) - x(1)) |
192 |
end if |
193 |
|
194 |
fpn = cs%c(n-1) + cs%b(n-1)*(cs%c(n-1) - cs%c(n-2))/(cs%b(n-2) + cs%b(n-1)) |
195 |
if (n > 3) then |
196 |
fpn = fpn + cs%b(n-1)*(cs%c(n-1) - cs%c(n-2) - (cs%b(n-2) + cs%b(n-1))* & |
197 |
(cs%c(n-2) - cs%c(n-3))/(cs%b(n-2) + cs%b(n-3)))/(x(n) - x(n-3)) |
198 |
end if |
199 |
|
200 |
! Calculate the right hand side and store it in C. |
201 |
|
202 |
cs%c(n) = 3.0_dp * (fpn - cs%c(n-1)) |
203 |
do i = n-1,2,-1 |
204 |
cs%c(i) = 3.0_dp * (cs%c(i) - cs%c(i-1)) |
205 |
end do |
206 |
cs%c(1) = 3.0_dp * (cs%c(1) - fp1) |
207 |
|
208 |
! Solve the tridiagonal system. |
209 |
|
210 |
do k = 2, n |
211 |
p = cs%b(k-1) / cs%d(k-1) |
212 |
cs%d(k) = cs%d(k) - p*cs%b(k-1) |
213 |
cs%c(k) = cs%c(k) - p*cs%c(k-1) |
214 |
end do |
215 |
cs%c(n) = cs%c(n) / cs%d(n) |
216 |
do k = n-1, 1, -1 |
217 |
cs%c(k) = (cs%c(k) - cs%b(k) * cs%c(k+1)) / cs%d(k) |
218 |
end do |
219 |
|
220 |
! Calculate the coefficients defining the spline. |
221 |
|
222 |
cs%d(1:n-1) = diff(cs%c) / (3.0_dp * h) |
223 |
cs%b(1:n-1) = diff_y / h - h * (cs%c(1:n-1) + h * cs%d(1:n-1)) |
224 |
cs%b(n) = cs%b(n-1) + h(n-1) * (2.0_dp*cs%c(n-1) + h(n-1)*3.0_dp*cs%d(n-1)) |
225 |
|
226 |
if (isUniform) then |
227 |
cs%dx_i = 1.0_dp / (x(2) - x(1)) |
228 |
endif |
229 |
|
230 |
return |
231 |
|
232 |
contains |
233 |
|
234 |
function diff(v) |
235 |
! Auxiliary function to compute the forward difference |
236 |
! of data stored in a vector v. |
237 |
|
238 |
implicit none |
239 |
real (kind = dp), dimension(:), intent(in) :: v |
240 |
real (kind = dp), dimension(size(v)-1) :: diff |
241 |
|
242 |
integer :: n |
243 |
|
244 |
n = size(v) |
245 |
diff = v(2:n) - v(1:n-1) |
246 |
return |
247 |
end function diff |
248 |
|
249 |
end subroutine newSpline |
250 |
|
251 |
subroutine deleteSpline(this) |
252 |
|
253 |
type(cubicSpline) :: this |
254 |
|
255 |
if(associated(this%d)) then |
256 |
deallocate(this%d) |
257 |
this%d => null() |
258 |
end if |
259 |
if(associated(this%c)) then |
260 |
deallocate(this%c) |
261 |
this%c => null() |
262 |
end if |
263 |
if(associated(this%b)) then |
264 |
deallocate(this%b) |
265 |
this%b => null() |
266 |
end if |
267 |
if(associated(this%y)) then |
268 |
deallocate(this%y) |
269 |
this%y => null() |
270 |
end if |
271 |
if(associated(this%x)) then |
272 |
deallocate(this%x) |
273 |
this%x => null() |
274 |
end if |
275 |
|
276 |
this%n = 0 |
277 |
|
278 |
end subroutine deleteSpline |
279 |
|
280 |
subroutine lookupNonuniformSpline(cs, xval, yval) |
281 |
|
282 |
implicit none |
283 |
|
284 |
type (cubicSpline), intent(in) :: cs |
285 |
real( kind = DP ), intent(in) :: xval |
286 |
real( kind = DP ), intent(out) :: yval |
287 |
real( kind = DP ) :: dx |
288 |
integer :: i, j |
289 |
! |
290 |
! Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains |
291 |
! or is nearest to xval. |
292 |
! |
293 |
j = cs%n - 1 |
294 |
|
295 |
do i = 0, cs%n - 2 |
296 |
|
297 |
if ( xval < cs%x(i+1) ) then |
298 |
j = i |
299 |
exit |
300 |
end if |
301 |
|
302 |
end do |
303 |
! |
304 |
! Evaluate the cubic polynomial. |
305 |
! |
306 |
dx = xval - cs%x(j) |
307 |
yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j))) |
308 |
|
309 |
return |
310 |
end subroutine lookupNonuniformSpline |
311 |
|
312 |
subroutine lookupUniformSpline(cs, xval, yval) |
313 |
|
314 |
implicit none |
315 |
|
316 |
type (cubicSpline), intent(in) :: cs |
317 |
real( kind = DP ), intent(in) :: xval |
318 |
real( kind = DP ), intent(out) :: yval |
319 |
real( kind = DP ) :: dx |
320 |
integer :: i, j |
321 |
! |
322 |
! Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains |
323 |
! or is nearest to xval. |
324 |
|
325 |
j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1)) |
326 |
|
327 |
dx = xval - cs%x(j) |
328 |
yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j))) |
329 |
|
330 |
return |
331 |
end subroutine lookupUniformSpline |
332 |
|
333 |
subroutine lookupUniformSpline1d(cs, xval, yval, dydx) |
334 |
|
335 |
implicit none |
336 |
|
337 |
type (cubicSpline), intent(in) :: cs |
338 |
real( kind = DP ), intent(in) :: xval |
339 |
real( kind = DP ), intent(out) :: yval, dydx |
340 |
real( kind = DP ) :: dx |
341 |
integer :: i, j |
342 |
|
343 |
! Find the interval J = [ cs%x(J), cs%x(J+1) ] that contains |
344 |
! or is nearest to xval. |
345 |
|
346 |
|
347 |
j = MAX(1, MIN(cs%n-1, int((xval-cs%x(1)) * cs%dx_i) + 1)) |
348 |
|
349 |
dx = xval - cs%x(j) |
350 |
yval = cs%y(j) + dx*(cs%b(j) + dx*(cs%c(j) + dx*cs%d(j))) |
351 |
|
352 |
dydx = cs%b(j) + dx*(2.0_dp * cs%c(j) + 3.0_dp * dx * cs%d(j)) |
353 |
|
354 |
return |
355 |
end subroutine lookupUniformSpline1d |
356 |
|
357 |
subroutine lookupSpline(cs, xval, yval) |
358 |
|
359 |
type (cubicSpline), intent(in) :: cs |
360 |
real( kind = DP ), intent(inout) :: xval |
361 |
real( kind = DP ), intent(inout) :: yval |
362 |
|
363 |
if (cs%isUniform) then |
364 |
call lookupUniformSpline(cs, xval, yval) |
365 |
else |
366 |
call lookupNonuniformSpline(cs, xval, yval) |
367 |
endif |
368 |
|
369 |
return |
370 |
end subroutine lookupSpline |
371 |
|
372 |
end module interpolation |