FV3 Bundle
cloud_interpolator.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 ! nf95 -r8 -g -I ~/regression/ia64/23-Jun-2005/CM2.1U_Control-1990_E1.k32pe/include/ -D_TEST_CLOUD_INTERPOLATOR -D_F95 cloud_interpolator.F90
20 
21 #define _FLATTEN(A) reshape((A), (/size((A))/) )
22 
24 #include <fms_platform.h>
25  implicit none
26  private
27 
28  public :: cld_ntrp_linear_cell_interp, cld_ntrp_locate_cell, cld_ntrp_get_cell_values
29 #ifdef _TEST_CLOUD_INTERPOLATOR
30  public :: cld_ntrp_expand_index, cld_ntrp_contract_indices
31 #endif
32 
33 ! Include variable "version" to be written to log file.
34 #include<file_version.h>
35 real, parameter :: tol = 10.0*epsilon(1.)
36 
37 CONTAINS
38 
39 !...............................................................................
40  _pure subroutine cld_ntrp_expand_index(Ic, ie, ier)
41  integer, intent(in) :: ic ! contacted index
42  integer, intent(out) :: ie(:) ! expanded list of indices
43  integer, intent(out) :: ier ! error flag (0=ok)
44 
45  integer j, nd
46 
47  ier = 0
48  nd = size(ie) ! dimension
49 
50  if(ic >= 2**nd) then
51  ie = -1
52  ier = 1 ! error
53  return
54  endif
55 
56  do j = 1, nd
57  ie(j) = mod(ic/2**(j-1), 2)
58  end do
59 
60  end subroutine cld_ntrp_expand_index
61 
62 !...............................................................................
63 !...............................................................................
64  _pure subroutine cld_ntrp_contract_indices(ie, Ic, ier)
65  integer, intent(in) :: ie(:) ! expanded list of indices
66  integer, intent(out) :: ic ! contacted index
67  integer, intent(out) :: ier ! error flag (0=ok)
68 
69  integer j, nd
70 
71  ier = 0
72  nd = size(ie) ! dimension
73 
74  ic = ie(nd)
75  do j = nd-1, 1, -1
76  ic = ic * 2
77  ic = ic + ie(j)
78  end do
79 
80  if(ic >= 2**nd) ier = 1
81 
82  end subroutine cld_ntrp_contract_indices
83 
84 
85 !...............................................................................
86 !...............................................................................
87  _pure subroutine cld_ntrp_linear_cell_interp(fvals, ts, f, ier)
88  real, intent(in) :: fvals(0:) ! values at the cell nodes
89  real, intent(in) :: ts(:) ! normalized [0,1]^nd cell coordinates
90  real, intent(out):: f ! interpolated value
91  integer, intent(out) :: ier ! error flag (0=ok)
92 
93  integer j, nd, ic, iflag
94  integer ie(size(fvals))
95  real basis
96 
97  ier = 0
98  f = 0.
99  nd = size(ts)
100  if(size(fvals) /= 2**nd) then
101  ier = 1
102  return
103  endif
104 
105  do ic = 0, 2**nd - 1
106  basis = 1.
107  call cld_ntrp_expand_index(ic, ie, iflag)
108  do j = 1, nd
109  basis = basis * ( (1.0-real(ie(j)))*(1.0-ts(j)) + real(ie(j))*ts(j) )
110  end do
111  f = f + fvals(ic)*basis
112  end do
113 
114  end subroutine cld_ntrp_linear_cell_interp
115 
116 !...............................................................................
117 !...............................................................................
118  _pure subroutine cld_ntrp_locate_cell(axis, x, index, ier)
119  real, intent(in) :: axis(:) ! axis
120  real, intent(in) :: x ! abscissae
121  integer, intent(out) :: index ! lower-left corner index
122  integer, intent(out) :: ier ! error flag (0=ok)
123 
124  logical down
125  integer n, index1, is
126  real axis_1, axis_n, axis_min, axis_max
127  ier = 0
128  index = -1
129  down = .false.
130  n = size(axis)
131  if(n < 2) then
132  ier = 3
133  return
134  endif
135  axis_1 = axis(1)
136  axis_n = axis(n)
137  axis_min = axis_1
138  axis_max = axis_n
139  if(axis_1 > axis_n) then
140  down = .true.
141  axis_min = axis_n
142  axis_max = axis_1
143  endif
144 
145  if(x < axis_min-tol) then
146  ier = 1
147  return
148  endif
149  if(x > axis_max+tol) then
150  ier = 2
151  return
152  endif
153 
154  index = floor(real(n-1)*(x - axis_1)/(axis_n-axis_1)) + 1
155  index = min(n-1, index)
156  index1 = index+1
157 
158  if(.NOT. down) then
159  if(axis(index) <= x+tol) then
160  if(x <= axis(index1)+tol) then
161  ! axis is uniform, or nearly so. Done!
162  return
163  else
164  ! increase index
165  is = index+1
166  do index = is, n-1
167  index1 = index+1
168  if(axis(index1) >= x-tol) return
169  enddo
170  endif
171  else
172  ! decrease index
173  is = index - 1
174  do index = is, 1, -1
175  if(axis(index) <= x+tol) return
176  enddo
177  endif
178  else
179  ! axis is pointing down
180  if(axis(index) >= x-tol) then
181  if(x >= axis(index1)-tol) then
182  ! axis is uniform, or nearly so. Done!
183  return
184  else
185  ! increase index
186  is = index + 1
187  do index = is, n-1
188  index1 = index+1
189  if(axis(index1) <= x+tol) return
190  enddo
191  endif
192  else
193  ! decrease index
194  is = index - 1
195  do index = is, 1, -1
196  if(axis(index) >= x-tol) return
197  enddo
198  endif
199  endif
200 
201  end subroutine cld_ntrp_locate_cell
202 
203 !...............................................................................
204 !...............................................................................
205  _pure subroutine cld_ntrp_get_flat_index(nsizes, indices, flat_index, ier)
206  integer, intent(in) :: nsizes(:) ! size of array along each axis
207  integer, intent(in) :: indices(:) ! cell indices
208  integer, intent(out) :: flat_index ! index into flattened array
209  integer, intent(out) :: ier ! error flag (0=ok)
210 
211  integer nd, id
212 
213  ier = 0
214  flat_index = -1
215  nd = size(nsizes)
216  if(nd /= size(indices)) then
217  ! size mismatch
218  ier = 1
219  return
220  endif
221 
222  flat_index = indices(nd)-1
223  do id = nd-1, 1, -1
224  flat_index = flat_index*nsizes(id) + indices(id)-1
225  enddo
226  flat_index = flat_index + 1
227 
228  end subroutine cld_ntrp_get_flat_index
229 
230 !...............................................................................
231 !...............................................................................
232  _pure subroutine cld_ntrp_get_cell_values(nsizes, fnodes, indices, fvals, ier)
233  integer, intent(in) :: nsizes(:) ! size of fnodes along each axis
234  real, intent(in) :: fnodes(:) ! flattened array of node values
235  integer, intent(in) :: indices(:) ! cell indices
236  real, intent(out) :: fvals(0:) ! returned array values in the cell
237  integer, intent(out) :: ier ! error flag (0=ok)
238 
239  integer id, nt, nd, flat_index, ic, iflag
240  integer, dimension(size(nsizes)) :: cell_indices, node_indices
241  ier = 0
242  fvals = 0.
243 
244  nd = size(nsizes)
245  if(nd /= size(indices)) then
246  ! size mismatch
247  ier = 1
248  return
249  endif
250  if(2**nd > size(fvals)) then
251  ! not enough elements to hold result
252  ier = 2
253  return
254  endif
255  nt = 1
256  do id = 1, nd
257  nt = nt * nsizes(id)
258  enddo
259  if(nt /= size(fnodes)) then
260  ! not enough node values
261  ier = 3
262  return
263  endif
264 
265  do ic = 0, 2**nd-1
266  call cld_ntrp_expand_index(ic, cell_indices, iflag)
267  node_indices = indices + cell_indices
268  call cld_ntrp_get_flat_index(nsizes, node_indices, flat_index, iflag)
269  fvals(ic) = fnodes(flat_index)
270  enddo
271 
272  end subroutine cld_ntrp_get_cell_values
273 
274 end MODULE cloud_interpolator_mod
275 !===============================================================================
276 
277 #ifdef _TEST_CLOUD_INTERPOLATOR
278 program test
280  implicit none
281 
282  call test_expansion_contraction
283  call test_linear_cell_interpolation
284  call test_cell_search
285  call test_get_node_values
286 
287  contains
288  subroutine test_expansion_contraction
289  integer ie1(4), ie2(4), ic, ier, idiff, j
290  ie1 = (/1,0,1,1/)
291  call cld_ntrp_contract_indices(ie1, ic, ier)
292  if(ier/=0) print *,'ERROR flag ier=', ier
293  call cld_ntrp_expand_index(ic, ie2, ier)
294  if(ier/=0) print *,'ERROR flag ier=', ier
295  idiff = 0
296  do j = 1, size(ie1)
297  idiff = idiff + abs(ie1(j)-ie2(j))
298  end do
299  if(idiff/=0) then
300  print *,'ERROR: contraction/expansion test failed (ie1/=ie2)'
301  endif
302  print *,'ie1 = ', ie1
303  print *,'ie2 = ', ie2
304  print *,'Ic = ', ic
305 
306  end subroutine test_expansion_contraction
307 
308  subroutine test_linear_cell_interpolation
309  integer, parameter :: nd = 3
310  real :: fvals(2**nd), ts(nd)
311  real :: fi, fx
312  integer ier
313  ! f = 1 + x + 2*y + 3*z
314  fvals = (/ 1., 2., 3., 4., 4., 5., 6., 7. /)
315  ts = (/ 0.1, 0.2, 0.3 /)
316  fx = 1. + ts(1) + 2*ts(2) + 3*ts(3)
317  call cld_ntrp_linear_cell_interp(fvals, ts, fi, ier)
318  if(ier/=0) print *,'ERROR flag ier=', ier
319  print *,'fi, fx = ', fi, fx
320  end subroutine test_linear_cell_interpolation
321 
322  subroutine test_cell_search
323  integer index, ier
324  integer, parameter :: n = 5
325  real :: axis1(n) = (/0., 0.1, 0.2, 0.3, 0.4/)
326  real :: axis2(n) = (/0., 0.01, 0.02, 0.03, 0.4/)
327  real :: axis3(n) = (/0.4, 0.3, 0.2, 0.1, 0./)
328  real :: axis4(n) = (/0.4, 0.03, 0.02, 0.01, 0./)
329  real x
330  integer :: ier_tot = 0
331 
332  print *,'axis1=', axis1
333  x = -0.0001
334  call cld_ntrp_locate_cell(axis1, x, index, ier)
335  print *,' x=',x, ' index=', index, ' ==? ', -1
336  ier_tot = ier_tot + abs(index - (-1))
337  x = 0.
338  call cld_ntrp_locate_cell(axis1, x, index, ier)
339  print *,' x=',x, ' index=', index, ' ==? ', 1
340  ier_tot = ier_tot + abs(index - (1))
341  x = 0.1
342  call cld_ntrp_locate_cell(axis1, x, index, ier)
343  print *, ' x=',x, ' index=', index, ' ==? ', 2
344  ier_tot = ier_tot + abs(index - (2))
345  x = 0.4
346  call cld_ntrp_locate_cell(axis1, x, index, ier)
347  print *,' x=',x, ' index=', index, ' ==? ', 4
348  ier_tot = ier_tot + abs(index - (4))
349  x = 0.40001
350  call cld_ntrp_locate_cell(axis1, x, index, ier)
351  print *,' x=',x, ' index=', index, ' ==? ', -1
352  ier_tot = ier_tot + abs(index - (-1))
353 
354  print *,'axis2=', axis1
355  x = -0.0001
356  call cld_ntrp_locate_cell(axis2, x, index, ier)
357  print *,' x=',x, ' index=', index, ' ==? ', -1
358  ier_tot = ier_tot + abs(index - (-1))
359  x = 0.
360  call cld_ntrp_locate_cell(axis2, x, index, ier)
361  print *,' x=',x, ' index=', index, ' ==? ', 1
362  ier_tot = ier_tot + abs(index - (1))
363  x = 0.1
364  call cld_ntrp_locate_cell(axis2, x, index, ier)
365  print *, ' x=',x, ' index=', index, ' ==? ', 4
366  ier_tot = ier_tot + abs(index - (4))
367  x = 0.4
368  call cld_ntrp_locate_cell(axis2, x, index, ier)
369  print *,' x=',x, ' index=', index, ' ==? ', 4
370  ier_tot = ier_tot + abs(index - (4))
371  x = 0.40001
372  call cld_ntrp_locate_cell(axis2, x, index, ier)
373  print *,' x=',x, ' index=', index, ' ==? ', -1
374  ier_tot = ier_tot + abs(index - (-1))
375 
376  print *,'axis3=', axis1
377  x = -0.0001
378  call cld_ntrp_locate_cell(axis3, x, index, ier)
379  print *,' x=',x, ' index=', index, ' ==? ', -1
380  ier_tot = ier_tot + abs(index - (-1))
381  x = 0.
382  call cld_ntrp_locate_cell(axis3, x, index, ier)
383  print *,' x=',x, ' index=', index, ' ==? ', 4
384  ier_tot = ier_tot + abs(index - (4))
385  x = 0.1
386  call cld_ntrp_locate_cell(axis3, x, index, ier)
387  print *, ' x=',x, ' index=', index, ' ==? ', 4
388  ier_tot = ier_tot + abs(index - (4))
389  x = 0.4
390  call cld_ntrp_locate_cell(axis3, x, index, ier)
391  print *,' x=',x, ' index=', index, ' ==? ', 1
392  ier_tot = ier_tot + abs(index - (1))
393  x = 0.40001
394  call cld_ntrp_locate_cell(axis3, x, index, ier)
395  print *,' x=',x, ' index=', index, ' ==? ', -1
396  ier_tot = ier_tot + abs(index - (-1))
397 
398  print *,'axis4=', axis1
399  x = -0.0001
400  call cld_ntrp_locate_cell(axis4, x, index, ier)
401  print *,' x=',x, ' index=', index, ' ==? ', -1
402  ier_tot = ier_tot + abs(index - (-1))
403  x = 0.
404  call cld_ntrp_locate_cell(axis4, x, index, ier)
405  print *,' x=',x, ' index=', index, ' ==? ', 4
406  ier_tot = ier_tot + abs(index - (4))
407  x = 0.1
408  call cld_ntrp_locate_cell(axis4, x, index, ier)
409  print *, ' x=',x, ' index=', index, ' ==? ', 1
410  ier_tot = ier_tot + abs(index - (1))
411  x = 0.4
412  call cld_ntrp_locate_cell(axis4, x, index, ier)
413  print *,' x=',x, ' index=', index, ' ==? ', 1
414  ier_tot = ier_tot + abs(index - (1))
415  x = 0.40001
416  call cld_ntrp_locate_cell(axis4, x, index, ier)
417  print *,' x=',x, ' index=', index, ' ==? ', -1
418  ier_tot = ier_tot + abs(index - (-1))
419 
420  print *,'Total error in test_cell_search: ', ier_tot
421 
422  end subroutine test_cell_search
423 
424  subroutine test_get_node_values
425  integer, parameter :: nd = 3, n1=6, n2=5, n3=4
426  real, dimension(n1, n2, n3) :: fnodes
427  real :: fvals(2**nd), fexact(2**nd)
428  real x, y, z
429  integer i, j, k, ier, indices(nd)
430  real :: error_tot = 0.
431  do k = 1, n3
432  do j = 1, n2
433  do i = 1, n1
434  x = 1* real(i-1)/real(n1-1)
435  y = 2* real(j-1)/real(n2-1)
436  z = 3* real(k-1)/real(n3-1)
437  fnodes(i,j,k) = x + y*z**2
438  enddo
439  enddo
440  enddo
441 
442  indices = (/1,1,1/)
443  call cld_ntrp_get_cell_values((/n1,n2,n3/), _flatten(fnodes), indices, fvals, ier)
444  fexact = (/0.0, 0.2, 0.0, 0.2, 0.0, 0.2, 0.5, 0.7/)
445  if(ier/=0) print *,'ERROR flag ier=', ier
446  print *,'indices ', indices
447  print *,'fvals=', fvals, ' ==? ', fexact
448  error_tot = error_tot + abs(sum(fvals - fexact))
449 
450  indices = (/5,4,2/)
451  call cld_ntrp_get_cell_values((/n1,n2,n3/), _flatten(fnodes), indices, fvals, ier)
452  fexact = (/2.3, 2.5, 2.8, 3.0, 6.8, 7.0, 8.8, 9.0/)
453  if(ier/=0) print *,'ERROR flag ier=', ier
454  print *,'indices ', indices
455  print *,'fvals=', fvals, ' ==? ', fexact
456  error_tot = error_tot + abs(sum(fvals - fexact))
457 
458  print *,'Total error in test_get_node_values: ', error_tot
459 
460 
461  end subroutine test_get_node_values
462 
463  end program test
464 
465 #endif
************************************************************************GNU Lesser General Public License **This file is part of the GFDL Flexible Modeling System(FMS). ! *! *FMS is free software 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 see< http:! ***********************************************************************subroutine READ_RECORD_CORE_(unit, field, nwords, data, start, axsiz) integer, intent(in) ::unit type(fieldtype), intent(in) ::field integer, intent(in) ::nwords MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in) ::start(:), axsiz(:) integer(SHORT_KIND) ::i2vals(nwords)!rab used in conjunction with transfer intrinsic to determine size of a variable integer(KIND=1) ::one_byte(8) integer ::word_sz!#ifdef __sgi integer(INT_KIND) ::ivals(nwords) real(FLOAT_KIND) ::rvals(nwords)!#else! integer ::ivals(nwords)! real ::rvals(nwords)!#endif real(DOUBLE_KIND) ::r8vals(nwords) pointer(ptr1, i2vals) pointer(ptr2, ivals) pointer(ptr3, rvals) pointer(ptr4, r8vals) if(mpp_io_stack_size< nwords) call mpp_io_set_stack_size(nwords) call mpp_error(FATAL, 'MPP_READ currently requires use_netCDF option') end subroutine READ_RECORD_CORE_ subroutine READ_RECORD_(unit, field, nwords, data, time_level, domain, position, tile_count, start_in, axsiz_in)!routine that is finally called by all mpp_read routines to perform the read!a non-netCDF record contains:! field ID! a set of 4 coordinates(is:ie, js:je) giving the data subdomain! a timelevel and a timestamp(=NULLTIME if field is static)! 3D real data(stored as 1D)!if you are using direct access I/O, the RECL argument to OPEN must be large enough for the above!in a global direct access file, record position on PE is given by %record.!Treatment of timestamp:! We assume that static fields have been passed without a timestamp.! Here that is converted into a timestamp of NULLTIME.! For non-netCDF fields, field is treated no differently, but is written! with a timestamp of NULLTIME. There is no check in the code to prevent! the user from repeatedly writing a static field. integer, intent(in) ::unit, nwords type(fieldtype), intent(in) ::field MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in), optional ::time_level type(domain2D), intent(in), optional ::domain integer, intent(in), optional ::position, tile_count integer, intent(in), optional ::start_in(:), axsiz_in(:) integer, dimension(size(field%axes(:))) ::start, axsiz integer ::tlevel !, subdomain(4) integer ::i, error, is, ie, js, je, isg, ieg, jsg, jeg type(domain2d), pointer ::io_domain=> tlevel if(PRESENT(start_in) .AND. PRESENT(axsiz_in)) then if(size(start(! the data domain and compute domain must refer to the subdomain being passed ! In this ! since that attempts to gather all data on PE size(field%axes(:)) axsiz(i)
#define min(a, b)
Definition: mosaic_util.h:32