-
Notifications
You must be signed in to change notification settings - Fork 140
/
Copy pathcloud_interpolator.F90
292 lines (252 loc) · 8.44 KB
/
cloud_interpolator.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
!***********************************************************************
!* GNU Lesser General Public License
!*
!* This file is part of the GFDL Flexible Modeling System (FMS).
!*
!* FMS is free software: you can redistribute it and/or modify it under
!* the terms of the GNU Lesser General Public License as published by
!* the Free Software Foundation, either version 3 of the License, or (at
!* your option) any later version.
!*
!* FMS is distributed in the hope that it will be useful, but WITHOUT
!* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
!* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
!* for more details.
!*
!* You should have received a copy of the GNU Lesser General Public
!* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
!***********************************************************************
#define _FLATTEN(A) reshape((A), (/size((A))/) )
!> @defgroup cloud_interpolator_mod cloud_interpolator_mod
!> @ingroup drifters
!! @brief Cloud interpolation routines for use in @ref drifters_mod
!> @addtogroup cloud_interpolator_mod
!> @{
MODULE cloud_interpolator_mod
#ifdef use_drifters
implicit none
private
public :: cld_ntrp_linear_cell_interp, cld_ntrp_locate_cell, cld_ntrp_get_cell_values
public :: cld_ntrp_expand_index, cld_ntrp_contract_indices
! Include variable "version" to be written to log file.
#include<file_version.h>
real, parameter :: tol = 10.0*epsilon(1.)
CONTAINS
!...............................................................................
!> Get expanded list of indices from contracted index
!> @param Ic contracted index
!> @param[out] ie(:) expanded list of indices
!> @param[out] ier error flag, non zero if operation unsuccessful
pure subroutine cld_ntrp_expand_index(Ic, ie, ier)
integer, intent(in) :: Ic
integer, intent(out) :: ie(:)
integer, intent(out) :: ier
integer j, nd
ier = 0
nd = size(ie) ! dimension
if(Ic >= 2**nd) then
ie = -1
ier = 1 ! error
return
endif
do j = 1, nd
ie(j) = mod(Ic/2**(j-1), 2)
end do
end subroutine cld_ntrp_expand_index
!...............................................................................
!> Contract list of indices to an single integer
!> @param ie(:) expanded list of indices
!> @param[out] Ic contracted index
!> @param[out] ier error flag, non zero if operation unsuccessful
pure subroutine cld_ntrp_contract_indices(ie, Ic, ier)
integer, intent(in) :: ie(:)
integer, intent(out) :: Ic
integer, intent(out) :: ier
integer j, nd
ier = 0
nd = size(ie) ! dimension
Ic = ie(nd)
do j = nd-1, 1, -1
Ic = Ic * 2
Ic = Ic + ie(j)
end do
if(Ic >= 2**nd) ier = 1
end subroutine cld_ntrp_contract_indices
!...............................................................................
!...............................................................................
!> Cloud interpolation for linear cells
!> @param fvals values at the cell nodes
!> @param ts normalized [0,1]^nd cell coordinates
!> @param[out] interpolated value
!> @param[out] error flag, non zero if unsucessful
pure subroutine cld_ntrp_linear_cell_interp(fvals, ts, f, ier)
real, intent(in) :: fvals(0:)
real, intent(in) :: ts(:)
real, intent(out):: f
integer, intent(out) :: ier
integer j, nd, Ic, iflag
integer ie(size(fvals))
real basis
ier = 0
f = 0.
nd = size(ts)
if(size(fvals) /= 2**nd) then
ier = 1
return
endif
do Ic = 0, 2**nd - 1
basis = 1.
call cld_ntrp_expand_index(Ic, ie, iflag)
do j = 1, nd
basis = basis * ( (1.0-real(ie(j)))*(1.0-ts(j)) + real(ie(j))*ts(j) )
end do
f = f + fvals(Ic)*basis
end do
end subroutine cld_ntrp_linear_cell_interp
!...............................................................................
!...............................................................................
pure subroutine cld_ntrp_locate_cell(axis, x, index, ier)
real, intent(in) :: axis(:) !< axis
real, intent(in) :: x !< abscissae
integer, intent(out) :: index !< lower-left corner index
integer, intent(out) :: ier !< error flag (0=ok)
logical down
integer n, index1, is
real axis_1, axis_n, axis_min, axis_max
ier = 0
index = -1
down = .FALSE.
n = size(axis)
if(n < 2) then
ier = 3
return
endif
axis_1 = axis(1)
axis_n = axis(n)
axis_min = axis_1
axis_max = axis_n
if(axis_1 > axis_n) then
down = .TRUE.
axis_min = axis_n
axis_max = axis_1
endif
if(x < axis_min-tol) then
ier = 1
return
endif
if(x > axis_max+tol) then
ier = 2
return
endif
index = floor(real(n-1)*(x - axis_1)/(axis_n-axis_1)) + 1
index = min(n-1, index)
index1 = index+1
if(.NOT. down) then
if(axis(index) <= x+tol) then
if(x <= axis(index1)+tol) then
! axis is uniform, or nearly so. Done!
return
else
! increase index
is = index+1
do index = is, n-1
index1 = index+1
if(axis(index1) >= x-tol) return
enddo
endif
else
! decrease index
is = index - 1
do index = is, 1, -1
if(axis(index) <= x+tol) return
enddo
endif
else
! axis is pointing down
if(axis(index) >= x-tol) then
if(x >= axis(index1)-tol) then
! axis is uniform, or nearly so. Done!
return
else
! increase index
is = index + 1
do index = is, n-1
index1 = index+1
if(axis(index1) <= x+tol) return
enddo
endif
else
! decrease index
is = index - 1
do index = is, 1, -1
if(axis(index) >= x-tol) return
enddo
endif
endif
end subroutine cld_ntrp_locate_cell
!...............................................................................
!...............................................................................
pure subroutine cld_ntrp_get_flat_index(nsizes, indices, flat_index, ier)
integer, intent(in) :: nsizes(:) !< size of array along each axis
integer, intent(in) :: indices(:) !< cell indices
integer, intent(out) :: flat_index !< index into flattened array
integer, intent(out) :: ier !< error flag (0=ok)
integer nd, id
ier = 0
flat_index = -1
nd = size(nsizes)
if(nd /= size(indices)) then
! size mismatch
ier = 1
return
endif
flat_index = indices(nd)-1
do id = nd-1, 1, -1
flat_index = flat_index*nsizes(id) + indices(id)-1
enddo
flat_index = flat_index + 1
end subroutine cld_ntrp_get_flat_index
!...............................................................................
!...............................................................................
pure subroutine cld_ntrp_get_cell_values(nsizes, fnodes, indices, fvals, ier)
integer, intent(in) :: nsizes(:) !< size of fnodes along each axis
real, intent(in) :: fnodes(:) !< flattened array of node values
integer, intent(in) :: indices(:) !< cell indices
real, intent(out) :: fvals(0:) !< returned array values in the cell
integer, intent(out) :: ier !< error flag (0=ok)
integer id, nt, nd, flat_index, Ic, iflag
integer, dimension(size(nsizes)) :: cell_indices, node_indices
ier = 0
fvals = 0.
nd = size(nsizes)
if(nd /= size(indices)) then
! size mismatch
ier = 1
return
endif
if(2**nd > size(fvals)) then
! not enough elements to hold result
ier = 2
return
endif
nt = 1
do id = 1, nd
nt = nt * nsizes(id)
enddo
if(nt /= size(fnodes)) then
! not enough node values
ier = 3
return
endif
do Ic = 0, 2**nd-1
call cld_ntrp_expand_index(Ic, cell_indices, iflag)
node_indices = indices + cell_indices
call cld_ntrp_get_flat_index(nsizes, node_indices, flat_index, iflag)
fvals(Ic) = fnodes(flat_index)
enddo
end subroutine cld_ntrp_get_cell_values
#endif
end MODULE cloud_interpolator_mod
!===============================================================================
!> @}
! close documentation grouping