1 |
chuckv |
930 |
!! |
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. Acknowledgement of the program authors must be made in any |
10 |
|
|
!! publication of scientific results based in part on use of the |
11 |
|
|
!! program. An acceptable form of acknowledgement is citation of |
12 |
|
|
!! the article in which the program was described (Matthew |
13 |
|
|
!! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
|
|
!! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
|
|
!! Parallel Simulation Engine for Molecular Dynamics," |
16 |
|
|
!! J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
|
|
!! |
18 |
|
|
!! 2. Redistributions of source code must retain the above copyright |
19 |
|
|
!! notice, this list of conditions and the following disclaimer. |
20 |
|
|
!! |
21 |
|
|
!! 3. Redistributions in binary form must reproduce the above copyright |
22 |
|
|
!! notice, this list of conditions and the following disclaimer in the |
23 |
|
|
!! documentation and/or other materials provided with the |
24 |
|
|
!! distribution. |
25 |
|
|
!! |
26 |
|
|
!! This software is provided "AS IS," without a warranty of any |
27 |
|
|
!! kind. All express or implied conditions, representations and |
28 |
|
|
!! warranties, including any implied warranty of merchantability, |
29 |
|
|
!! fitness for a particular purpose or non-infringement, are hereby |
30 |
|
|
!! excluded. The University of Notre Dame and its licensors shall not |
31 |
|
|
!! be liable for any damages suffered by licensee as a result of |
32 |
|
|
!! using, modifying or distributing the software or its |
33 |
|
|
!! derivatives. In no event will the University of Notre Dame or its |
34 |
|
|
!! licensors be liable for any lost revenue, profit or data, or for |
35 |
|
|
!! direct, indirect, special, consequential, incidental or punitive |
36 |
|
|
!! damages, however caused and regardless of the theory of liability, |
37 |
|
|
!! arising out of the use of or inability to use software, even if the |
38 |
|
|
!! University of Notre Dame has been advised of the possibility of |
39 |
|
|
!! such damages. |
40 |
|
|
!! |
41 |
|
|
!! |
42 |
|
|
!! interpolation.F90 |
43 |
|
|
!! |
44 |
|
|
!! Created by Charles F. Vardeman II on 03 Apr 2006. |
45 |
|
|
!! |
46 |
|
|
!! PURPOSE: Generic Spline interplelation routines. These routines assume that we are on a uniform grid for |
47 |
|
|
!! precomputation of spline parameters. |
48 |
|
|
!! |
49 |
|
|
!! @author Charles F. Vardeman II |
50 |
|
|
!! @version $Id: interpolation.F90,v 1.1 2006-04-12 21:15:17 chuckv Exp $ |
51 |
|
|
|
52 |
|
|
|
53 |
|
|
module INTERPOLATION |
54 |
|
|
use definitions |
55 |
|
|
use status |
56 |
|
|
implicit none |
57 |
|
|
PRIVATE |
58 |
|
|
|
59 |
|
|
character(len = statusMsgSize) :: errMSG |
60 |
|
|
|
61 |
|
|
type, public :: splineType |
62 |
|
|
private |
63 |
|
|
integer :: npoints = 0 |
64 |
|
|
real(kind=dp) :: delta_x |
65 |
|
|
real(kind=dp) :: range |
66 |
|
|
real(kind=dp) :: range_inv |
67 |
|
|
real (kind=dp), pointer,dimension(:) :: xa => null() |
68 |
|
|
real (kind=dp), pointer,dimension(:) :: ya => null() |
69 |
|
|
real (kind=dp), pointer,dimension(:) :: yppa => null() |
70 |
|
|
end type splineType |
71 |
|
|
|
72 |
|
|
type, public :: multiSplineType |
73 |
|
|
private |
74 |
|
|
integer :: npoints = 0 |
75 |
|
|
integer :: nfuncs = 0 |
76 |
|
|
|
77 |
|
|
integer :: npoints = 0 |
78 |
|
|
real(kind=dp) :: delta_x |
79 |
|
|
real(kind=dp) :: range |
80 |
|
|
real(kind=dp) :: range_inv |
81 |
|
|
real (kind=dp), pointer,dimension(:) :: xa => null() |
82 |
|
|
real (kind=dp), pointer,dimension(:,:) :: ya => null() |
83 |
|
|
real (kind=dp), pointer,dimension(:,:) :: yppa => null() |
84 |
|
|
end type splineType |
85 |
|
|
|
86 |
|
|
|
87 |
|
|
interface splineLookup |
88 |
|
|
module procedure multiSplint |
89 |
|
|
module procedure splintd |
90 |
|
|
module procedure splintd1 |
91 |
|
|
module procedure splintd2 |
92 |
|
|
end interface |
93 |
|
|
|
94 |
|
|
public :: splint |
95 |
|
|
public :: newSpline |
96 |
|
|
public :: newMultiSpline |
97 |
|
|
public :: deleteSpline |
98 |
|
|
public :: deleteMultiSpline |
99 |
|
|
|
100 |
|
|
|
101 |
|
|
contains |
102 |
|
|
|
103 |
|
|
!! mySpline is pointer to spline type, nx number of data points, |
104 |
|
|
!! xa tabulated x values and ya respective values for xa, yp1 |
105 |
|
|
!! is value for derivative at first point and ypn is value |
106 |
|
|
!! for derivative at point n. |
107 |
|
|
subroutine newSpline(thisSpline,nx, xa, ya, yp1, ypn, boundary) |
108 |
|
|
|
109 |
|
|
! yp1 and ypn are the first derivatives of y at the two endpoints |
110 |
|
|
! if boundary is 'L' the lower derivative is used |
111 |
|
|
! if boundary is 'U' the upper derivative is used |
112 |
|
|
! if boundary is 'B' then both derivatives are used |
113 |
|
|
! if boundary is anything else, then both derivatives are assumed to be 0 |
114 |
|
|
|
115 |
|
|
|
116 |
|
|
type (splineType), intent(inout) :: thisSpline |
117 |
|
|
|
118 |
|
|
|
119 |
|
|
real( kind = DP ), pointer, dimension(:) :: xa |
120 |
|
|
real( kind = DP ), pointer, dimension(:) :: ya |
121 |
|
|
real( kind = DP ), dimension(size(xa)) :: u |
122 |
|
|
real( kind = DP ) :: yp1,ypn,un,qn,sig,p |
123 |
|
|
character(len=*) :: boundary |
124 |
|
|
integer :: nx, i, k, max_array_size |
125 |
|
|
integer :: alloc_error |
126 |
|
|
|
127 |
|
|
alloc_error = 0 |
128 |
|
|
|
129 |
|
|
if (thisSpline%npoints/=0) then |
130 |
|
|
call handleWarning("INTERPOLATION:newSpline",& |
131 |
|
|
"Type has already been created") |
132 |
|
|
call deleteSpline(thisSpline) |
133 |
|
|
end if |
134 |
|
|
|
135 |
|
|
|
136 |
|
|
! make sure the sizes match |
137 |
|
|
if ((nx /= size(xa)) .or. (nx /= size(ya))) then |
138 |
|
|
call handleWarning("INTERPOLATION:newSpline","Array size mismatch") |
139 |
|
|
end if |
140 |
|
|
|
141 |
|
|
|
142 |
|
|
thisSpline%npoints = nx |
143 |
|
|
allocate(thisSpline%yppa(nx),stat=alloc_error) |
144 |
|
|
if(alloc_error .ne. 0) call handleWarning("INTERPOLATION:newSpline",& |
145 |
|
|
"Error in allocating storage for yppa") |
146 |
|
|
|
147 |
|
|
thisSpline%xa => xa |
148 |
|
|
thisSpline%ya => ya |
149 |
|
|
|
150 |
|
|
|
151 |
|
|
|
152 |
|
|
|
153 |
|
|
if ((boundary.eq.'l').or.(boundary.eq.'L').or. & |
154 |
|
|
(boundary.eq.'b').or.(boundary.eq.'B')) then |
155 |
|
|
thisSpline%yppa(1) = -0.5E0_DP |
156 |
|
|
u(1) = (3.0E0_DP/(xa(2)-xa(1)))*((ya(2)-& |
157 |
|
|
ya(1))/(xa(2)-xa(1))-yp1) |
158 |
|
|
else |
159 |
|
|
thisSpline%yppa(1) = 0.0E0_DP |
160 |
|
|
u(1) = 0.0E0_DP |
161 |
|
|
endif |
162 |
|
|
|
163 |
|
|
do i = 2, nx - 1 |
164 |
|
|
sig = (thisSpline%xa(i) - thisSpline%xa(i-1)) / (thisSpline%xa(i+1) - thisSpline%xa(i-1)) |
165 |
|
|
p = sig * thisSpline%yppa(i-1) + 2.0E0_DP |
166 |
|
|
thisSpline%yppa(i) = (sig - 1.0E0_DP) / p |
167 |
|
|
u(i) = (6.0E0_DP*((thisSpline%ya(i+1)-thisSpline%ya(i))/(thisSpline%xa(i+1)-thisSpline%xa(i)) - & |
168 |
|
|
(thisSpline%ya(i)-thisSpline%ya(i-1))/(thisSpline%xa(i)-thisSpline%xa(i-1)))/ & |
169 |
|
|
(thisSpline%xa(i+1)-thisSpline%xa(i-1)) - sig * u(i-1))/p |
170 |
|
|
enddo |
171 |
|
|
|
172 |
|
|
if ((boundary.eq.'u').or.(boundary.eq.'U').or. & |
173 |
|
|
(boundary.eq.'b').or.(boundary.eq.'B')) then |
174 |
|
|
qn = 0.5E0_DP |
175 |
|
|
un = (3.0E0_DP/(thisSpline%xa(nx)-thisSpline%xa(nx-1)))* & |
176 |
|
|
(ypn-(thisSpline%ya(nx)-thisSpline%ya(nx-1))/(thisSpline%xa(nx)-thisSpline%xa(nx-1))) |
177 |
|
|
else |
178 |
|
|
qn = 0.0E0_DP |
179 |
|
|
un = 0.0E0_DP |
180 |
|
|
endif |
181 |
|
|
|
182 |
|
|
thisSpline%yppa(nx)=(un-qn*u(nx-1))/(qn*thisSpline%yppa(nx-1)+1.0E0_DP) |
183 |
|
|
|
184 |
|
|
do k = nx-1, 1, -1 |
185 |
|
|
thisSpline%yppa(k)=thisSpline%yppa(k)*thisSpline%yppa(k+1)+u(k) |
186 |
|
|
enddo |
187 |
|
|
|
188 |
|
|
end subroutine newSpline |
189 |
|
|
|
190 |
|
|
subroutine deleteSpline(thisSpline) |
191 |
|
|
type(splineType) :: thisSpline |
192 |
|
|
|
193 |
|
|
|
194 |
|
|
|
195 |
|
|
if(associated(thisSpline%xa)) then |
196 |
|
|
deallocate(thisSpline%xa) |
197 |
|
|
thisSpline%xa => null() |
198 |
|
|
end if |
199 |
|
|
if(associated(thisSpline%ya)) then |
200 |
|
|
deallocate(thisSpline%ya) |
201 |
|
|
thisSpline%ya => null() |
202 |
|
|
end if |
203 |
|
|
if(associated(thisSpline%yppa)) then |
204 |
|
|
deallocate(thisSpline%yppa) |
205 |
|
|
thisSpline%yppa => null() |
206 |
|
|
end if |
207 |
|
|
|
208 |
|
|
thisSpline%npoints=0 |
209 |
|
|
|
210 |
|
|
end subroutine deleteSpline |
211 |
|
|
|
212 |
|
|
subroutine splintd2(thisSpline, x, y, dy, d2y) |
213 |
|
|
type(splineType) :: thisSpline |
214 |
|
|
real( kind = DP ), intent(in) :: x |
215 |
|
|
real( kind = DP ), intent(out) :: y,dy,d2y |
216 |
|
|
|
217 |
|
|
|
218 |
|
|
real( kind = DP ) :: del, h, a, b, c, d |
219 |
|
|
integer :: j |
220 |
|
|
|
221 |
|
|
! this spline code assumes that the x points are equally spaced |
222 |
|
|
! do not attempt to use this code if they are not. |
223 |
|
|
|
224 |
|
|
|
225 |
|
|
! find the closest point with a value below our own: |
226 |
|
|
j = FLOOR(real((thisSpline%npoints-1),kind=dp) * & |
227 |
|
|
(x - thisSpline%xa(1)) / (thisSpline%xa(thisSpline%npoints) - thisSpline%xa(1))) + 1 |
228 |
|
|
|
229 |
|
|
! check to make sure we're inside the spline range: |
230 |
|
|
if ((j.gt.thisSpline%npoints).or.(j.lt.1)) then |
231 |
|
|
write(errMSG,*) "EAM_splint: x is outside bounds of spline: ",x,j |
232 |
|
|
call handleError("INTERPOLATION::SPLINT2d",errMSG) |
233 |
|
|
endif |
234 |
|
|
! check to make sure we haven't screwed up the calculation of j: |
235 |
|
|
if ((x.lt.thisSpline%xa(j)).or.(x.gt.thisSpline%xa(j+1))) then |
236 |
|
|
if (j.ne.thisSpline%npoints) then |
237 |
|
|
write(errMSG,*) "EAM_splint:",x," x is outside bounding range" |
238 |
|
|
call handleError("INTERPOLATION::SPLINT2d",errMSG) |
239 |
|
|
endif |
240 |
|
|
endif |
241 |
|
|
|
242 |
|
|
del = thisSpline%xa(j+1) - x |
243 |
|
|
h = thisSpline%xa(j+1) - thisSpline%xa(j) |
244 |
|
|
|
245 |
|
|
a = del / h |
246 |
|
|
b = 1.0E0_DP - a |
247 |
|
|
c = a*(a*a - 1.0E0_DP)*h*h/6.0E0_DP |
248 |
|
|
d = b*(b*b - 1.0E0_DP)*h*h/6.0E0_DP |
249 |
|
|
|
250 |
|
|
y = a*thisSpline%ya(j) + b*thisSpline%ya(j+1) + c*thisSpline%yppa(j) + d*thisSpline%yppa(j+1) |
251 |
|
|
|
252 |
|
|
dy = (thisSpline%ya(j+1)-thisSpline%ya(j))/h & |
253 |
|
|
- (3.0E0_DP*a*a - 1.0E0_DP)*h*thisSpline%yppa(j)/6.0E0_DP & |
254 |
|
|
+ (3.0E0_DP*b*b - 1.0E0_DP)*h*thisSpline%yppa(j+1)/6.0E0_DP |
255 |
|
|
|
256 |
|
|
|
257 |
|
|
d2y = a*thisSpline%yppa(j) + b*thisSpline%yppa(j+1) |
258 |
|
|
|
259 |
|
|
|
260 |
|
|
end subroutine splintd2 |
261 |
|
|
subroutine splintd1(thisSpline, x, y, dy) |
262 |
|
|
type(splineType) :: thisSpline |
263 |
|
|
real( kind = DP ), intent(in) :: x |
264 |
|
|
real( kind = DP ), intent(out) :: y,dy |
265 |
|
|
|
266 |
|
|
|
267 |
|
|
real( kind = DP ) :: del, h, a, b, c, d |
268 |
|
|
integer :: j |
269 |
|
|
|
270 |
|
|
! this spline code assumes that the x points are equally spaced |
271 |
|
|
! do not attempt to use this code if they are not. |
272 |
|
|
|
273 |
|
|
|
274 |
|
|
! find the closest point with a value below our own: |
275 |
|
|
j = FLOOR(real((thisSpline%npoints-1),kind=dp) *& |
276 |
|
|
(x - thisSpline%xa(1)) / (thisSpline%xa(thisSpline%npoints) - thisSpline%xa(1))) + 1 |
277 |
|
|
|
278 |
|
|
! check to make sure we're inside the spline range: |
279 |
|
|
if ((j.gt.thisSpline%npoints).or.(j.lt.1)) then |
280 |
|
|
write(errMSG,*) "EAM_splint: x is outside bounds of spline: ",x,j |
281 |
|
|
call handleError("INTERPOLATION::SPLINT2d",errMSG) |
282 |
|
|
endif |
283 |
|
|
! check to make sure we haven't screwed up the calculation of j: |
284 |
|
|
if ((x.lt.thisSpline%xa(j)).or.(x.gt.thisSpline%xa(j+1))) then |
285 |
|
|
if (j.ne.thisSpline%npoints) then |
286 |
|
|
write(errMSG,*) "EAM_splint:",x," x is outside bounding range" |
287 |
|
|
call handleError("INTERPOLATION::SPLINT2d",errMSG) |
288 |
|
|
endif |
289 |
|
|
endif |
290 |
|
|
|
291 |
|
|
del = thisSpline%xa(j+1) - x |
292 |
|
|
h = thisSpline%xa(j+1) - thisSpline%xa(j) |
293 |
|
|
|
294 |
|
|
a = del / h |
295 |
|
|
b = 1.0E0_DP - a |
296 |
|
|
c = a*(a*a - 1.0E0_DP)*h*h/6.0E0_DP |
297 |
|
|
d = b*(b*b - 1.0E0_DP)*h*h/6.0E0_DP |
298 |
|
|
|
299 |
|
|
y = a*thisSpline%ya(j) + b*thisSpline%ya(j+1) + c*thisSpline%yppa(j) + d*thisSpline%yppa(j+1) |
300 |
|
|
|
301 |
|
|
dy = (thisSpline%ya(j+1)-thisSpline%ya(j))/h & |
302 |
|
|
- (3.0E0_DP*a*a - 1.0E0_DP)*h*thisSpline%yppa(j)/6.0E0_DP & |
303 |
|
|
+ (3.0E0_DP*b*b - 1.0E0_DP)*h*thisSpline%yppa(j+1)/6.0E0_DP |
304 |
|
|
|
305 |
|
|
|
306 |
|
|
|
307 |
|
|
|
308 |
|
|
|
309 |
|
|
end subroutine splintd1 |
310 |
|
|
subroutine splintd(thisSpline, x, y) |
311 |
|
|
type(splineType) :: thisSpline |
312 |
|
|
real( kind = DP ), intent(in) :: x |
313 |
|
|
real( kind = DP ), intent(out) :: y |
314 |
|
|
|
315 |
|
|
|
316 |
|
|
real( kind = DP ) :: del, h, a, b, c, d |
317 |
|
|
integer :: j |
318 |
|
|
|
319 |
|
|
! this spline code assumes that the x points are equally spaced |
320 |
|
|
! do not attempt to use this code if they are not. |
321 |
|
|
|
322 |
|
|
|
323 |
|
|
! find the closest point with a value below our own: |
324 |
|
|
j = FLOOR(real((thisSpline%npoints-1),kind=dp) * & |
325 |
|
|
(x - thisSpline%xa(1)) / (thisSpline%xa(thisSpline%npoints) - thisSpline%xa(1))) + 1 |
326 |
|
|
|
327 |
|
|
! check to make sure we're inside the spline range: |
328 |
|
|
if ((j.gt.thisSpline%npoints).or.(j.lt.1)) then |
329 |
|
|
write(errMSG,*) "EAM_splint: x is outside bounds of spline: ",x,j |
330 |
|
|
call handleError("INTERPOLATION::SPLINT2d",errMSG) |
331 |
|
|
endif |
332 |
|
|
! check to make sure we haven't screwed up the calculation of j: |
333 |
|
|
if ((x.lt.thisSpline%xa(j)).or.(x.gt.thisSpline%xa(j+1))) then |
334 |
|
|
if (j.ne.thisSpline%npoints) then |
335 |
|
|
write(errMSG,*) "EAM_splint:",x," x is outside bounding range" |
336 |
|
|
call handleError("INTERPOLATION::SPLINT2d",errMSG) |
337 |
|
|
endif |
338 |
|
|
endif |
339 |
|
|
|
340 |
|
|
del = thisSpline%xa(j+1) - x |
341 |
|
|
h = thisSpline%xa(j+1) - thisSpline%xa(j) |
342 |
|
|
|
343 |
|
|
a = del / h |
344 |
|
|
b = 1.0E0_DP - a |
345 |
|
|
c = a*(a*a - 1.0E0_DP)*h*h/6.0E0_DP |
346 |
|
|
d = b*(b*b - 1.0E0_DP)*h*h/6.0E0_DP |
347 |
|
|
|
348 |
|
|
y = a*thisSpline%ya(j) + b*thisSpline%ya(j+1) + c*thisSpline%yppa(j) + d*thisSpline%yppa(j+1) |
349 |
|
|
|
350 |
|
|
end subroutine splintd |
351 |
|
|
|
352 |
|
|
|
353 |
|
|
end module INTERPOLATION |