FV3 Bundle
axis_utils.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 
21  !
22  !<CONTACT EMAIL="Matthew.Harrison@noaa.gov">M.J. Harrison</CONTACT>
23  !
24  !<REVIEWER EMAIL="Bruce.Wyman@noaa.gov">Bruce Wyman</REVIEWER>
25  !
26 
27  !<OVERVIEW>
28  ! A set of utilities for manipulating axes and extracting axis
29  ! attributes
30  !</OVERVIEW>
31 
32  !<DESCRIPTION>
33  !
34  ! subroutine get_axis_cart(axis,cart) : Returns X,Y,Z or T cartesian attribute
35  ! subroutine get_axis_bounds(axis,axis_bound,axes) : Return axis_bound either from an array of
36  ! available axes, or defined based on axis mid-points
37  ! function get_axis_modulo : Returns true if axis has the modulo attribute
38  ! function get_axis_fold : Returns is axis is folded at a boundary (non-standard meta-data)
39  ! function lon_in_range : Returns lon_strt <= longitude <= lon_strt+360
40  ! subroutine tranlon : Returns monotonic array of longitudes s.t., lon_strt <= lon(:) <= lon_strt+360.
41  ! subroutine nearest_index : Return index of nearest point along axis
42  !
43  !</DESCRIPTION>
44  !
45 
47  mpp_get_atts, mpp_get_axis_data, mpp_modify_meta, &
48  mpp_get_att_name, mpp_get_att_type, mpp_get_att_char, &
49  mpp_get_att_length, mpp_get_axis_bounds
50  use mpp_mod, only: mpp_error, fatal, stdout
51  use fms_mod, only: lowercase, string_array_index, fms_error_handler
52 
53  implicit none
54 
55 # include <netcdf.inc>
56 
59 
60  private
61 
62  integer, parameter :: maxatts = 100
63  real, parameter :: epsln= 1.e-10
64  real, parameter :: fp5 = 0.5, f360 = 360.0
65 
66 ! Include variable "version" to be written to log file.
67 #include<file_version.h>
68 
69  interface interp_1d
70  module procedure interp_1d_1d
71  module procedure interp_1d_2d
72  module procedure interp_1d_3d
73  end interface
74 
75 contains
76 
77 
78  subroutine get_axis_cart(axis, cart)
79 
80  type(axistype), intent(in) :: axis
81  character(len=1), intent(out) :: cart
82  character(len=1) :: axis_cart
83  character(len=16), dimension(2) :: lon_names, lat_names
84  character(len=16), dimension(3) :: z_names
85  character(len=16), dimension(2) :: t_names
86  character(len=16), dimension(3) :: lon_units, lat_units
87  character(len=8) , dimension(4) :: z_units
88  character(len=3) , dimension(6) :: t_units
89  character(len=32) :: name
90  integer :: i,j
91 
92  lon_names = (/'lon','x '/)
93  lat_names = (/'lat','y '/)
94  z_names = (/'depth ','height','z '/)
95  t_names = (/'time','t '/)
96  lon_units = (/'degrees_e ', 'degrees_east', 'degreese '/)
97  lat_units = (/'degrees_n ', 'degrees_north', 'degreesn '/)
98  z_units = (/'cm ','m ','pa ','hpa'/)
99  t_units = (/'sec', 'min','hou','day','mon','yea'/)
100 
101  call mpp_get_atts(axis,cartesian=axis_cart)
102  cart = 'N'
103 
104  if ( lowercase(axis_cart) == 'x' ) cart = 'X'
105  if ( lowercase(axis_cart) == 'y' ) cart = 'Y'
106  if ( lowercase(axis_cart) == 'z' ) cart = 'Z'
107  if ( lowercase(axis_cart) == 't' ) cart = 'T'
108 
109  if (cart /= 'X' .and. cart /= 'Y' .and. cart /= 'Z' .and. cart /= 'T') then
110  call mpp_get_atts(axis,name=name)
111  name = lowercase(name)
112  do i=1,size(lon_names(:))
113  if (trim(name(1:3)) == trim(lon_names(i))) cart = 'X'
114  enddo
115  do i=1,size(lat_names(:))
116  if (trim(name(1:3)) == trim(lat_names(i))) cart = 'Y'
117  enddo
118  do i=1,size(z_names(:))
119  if (trim(name) == trim(z_names(i))) cart = 'Z'
120  enddo
121  do i=1,size(t_names(:))
122  if (trim(name) == t_names(i)) cart = 'T'
123  enddo
124  end if
125 
126  if (cart /= 'X' .and. cart /= 'Y' .and. cart /= 'Z' .and. cart /= 'T') then
127  call mpp_get_atts(axis,units=name)
128  name = lowercase(name)
129  do i=1,size(lon_units(:))
130  if (trim(name) == trim(lon_units(i))) cart = 'X'
131  enddo
132  do i=1,size(lat_units(:))
133  if (trim(name) == trim(lat_units(i))) cart = 'Y'
134  enddo
135  do i=1,size(z_units(:))
136  if (trim(name) == trim(z_units(i))) cart = 'Z'
137  enddo
138  do i=1,size(t_units(:))
139  if (name(1:3) == trim(t_units(i))) cart = 'T'
140  enddo
141  end if
142 
143  return
144 
145  end subroutine get_axis_cart
146 
147 
148  subroutine get_axis_bounds(axis,axis_bound,axes,bnd_name,err_msg)
150  type(axistype), intent(in) :: axis
151  type(axistype), intent(inout) :: axis_bound
152  type(axistype), intent(in), dimension(:) :: axes
153  character(len=*), intent(inout), optional :: bnd_name
154  character(len=*), intent(out), optional :: err_msg
155 
156  real, dimension(:), allocatable :: data, tmp
157 
158  integer :: i, len
159  character(len=128) :: name, units
160  character(len=256) :: longname
161  character(len=1) :: cartesian
162  logical :: bounds_found
163 
164  if(present(err_msg)) then
165  err_msg = ''
166  endif
167  axis_bound = default_axis
168  call mpp_get_atts(axis,units=units,longname=longname,&
169  cartesian=cartesian, len=len)
170  if(len .LE. 0) return
171  allocate(data(len+1))
172 
173  bounds_found = mpp_get_axis_bounds(axis, data, name=name)
174  longname = trim(longname)//' bounds'
175 
176  if(.not.bounds_found .and. len>1 ) then
177  ! The following calculation can not be done for len=1
178  call mpp_get_atts(axis,name=name)
179  name = trim(name)//'_bnds'
180  allocate(tmp(len))
181  call mpp_get_axis_data(axis,tmp)
182  do i=2,len
183  data(i)= tmp(i-1)+fp5*(tmp(i)-tmp(i-1))
184  enddo
185  data(1)= tmp(1)- fp5*(tmp(2)-tmp(1))
186  if (abs(data(1)) < epsln) data(1) = 0.0
187  data(len+1)= tmp(len)+ fp5*(tmp(len)-tmp(len-1))
188  if (data(1) == 0.0) then
189  if (abs(data(len+1)-360.) > epsln) data(len+1)=360.0
190  endif
191  endif
192  if(bounds_found .OR. len>1) then
193  call mpp_modify_meta(axis_bound,name=name,units=units,longname=&
194  longname,cartesian=cartesian,data=data)
195  endif
196  if(allocated(tmp)) deallocate(tmp)
197  deallocate(data)
198 
199  return
200  end subroutine get_axis_bounds
201 
202  function get_axis_modulo(axis)
204  type(axistype) :: axis
205  logical :: get_axis_modulo
206  integer :: natt, i
207  type(atttype), dimension(:), allocatable :: atts
208 
209 
210  call mpp_get_atts(axis,natts=natt)
211  allocate(atts(natt))
212  call mpp_get_atts(axis,atts=atts)
213 
214  get_axis_modulo=.false.
215  do i = 1,natt
216  if (lowercase(trim(mpp_get_att_name(atts(i)))) == 'modulo') get_axis_modulo = .true.
217  enddo
218 
219  deallocate(atts)
220 
221  return
222  end function get_axis_modulo
223 
224  function get_axis_modulo_times(axis, tbeg, tend)
226  logical :: get_axis_modulo_times
227  type(axistype), intent(in) :: axis
228  character(len=*), intent(out) :: tbeg, tend
229  integer :: natt, i
230  type(atttype), dimension(:), allocatable :: atts
231  logical :: found_tbeg, found_tend
232 
233  call mpp_get_atts(axis,natts=natt)
234  allocate(atts(natt))
235  call mpp_get_atts(axis,atts=atts)
236 
237  found_tbeg = .false.
238  found_tend = .false.
239 
240  do i = 1,natt
241  if(lowercase(trim(mpp_get_att_name(atts(i)))) == 'modulo_beg') then
242  if(mpp_get_att_length(atts(i)) > len(tbeg)) then
243  call mpp_error(fatal,'error in get: len(tbeg) too small to hold attribute')
244  endif
245  tbeg = trim(mpp_get_att_char(atts(i)))
246  found_tbeg = .true.
247  endif
248  if(lowercase(trim(mpp_get_att_name(atts(i)))) == 'modulo_end') then
249  if(mpp_get_att_length(atts(i)) > len(tend)) then
250  call mpp_error(fatal,'error in get: len(tend) too small to hold attribute')
251  endif
252  tend = trim(mpp_get_att_char(atts(i)))
253  found_tend = .true.
254  endif
255  enddo
256 
257  if(found_tbeg .and. .not.found_tend) then
258  call mpp_error(fatal,'error in get: Found modulo_beg but not modulo_end')
259  endif
260  if(.not.found_tbeg .and. found_tend) then
261  call mpp_error(fatal,'error in get: Found modulo_end but not modulo_beg')
262  endif
263 
264  get_axis_modulo_times = found_tbeg
265 
266  end function get_axis_modulo_times
267 
268  function get_axis_fold(axis)
270  type(axistype) :: axis
271  logical :: get_axis_fold
272  integer :: natt, i
273  type(atttype), dimension(:), allocatable :: atts
274 
275 
276  call mpp_get_atts(axis,natts=natt)
277  allocate(atts(natt))
278  call mpp_get_atts(axis,atts=atts)
279 
280  get_axis_fold=.false.
281  do i = 1,natt
282  if (mpp_get_att_char(atts(i)) == 'fold_top') get_axis_fold = .true.
283  enddo
284 
285  deallocate(atts)
286 
287  return
288  end function get_axis_fold
289 
290  function lon_in_range(lon, l_strt)
291  real :: lon, l_strt, lon_in_range, l_end
292 
293  lon_in_range = lon
294  l_end = l_strt+360.
295 
296  if (abs(lon_in_range - l_strt) < 1.e-4) then
297  lon_in_range = l_strt
298  return
299  endif
300 
301  if (abs(lon_in_range - l_end) < 1.e-4) then
302  lon_in_range = l_strt
303  return
304  endif
305 
306  do
307  if (lon_in_range < l_strt) then
309  else if (lon_in_range > l_end) then
311  else
312  exit
313  end if
314  end do
315 
316  end function lon_in_range
317 
318  subroutine tranlon(lon, lon_start, istrt)
320  ! returns array of longitudes s.t. lon_strt <= lon < lon_strt+360.
321  ! also, the first istrt-1 entries are moved to the end of the array
322  !
323  ! e.g.
324  ! lon = 0 1 2 3 4 5 ... 358 359; lon_strt = 3 ==>
325  ! tranlon = 3 4 5 6 7 8 ... 359 360 361 362; istrt = 4
326 
327  real, intent(inout), dimension(:) :: lon
328  real, intent(in) :: lon_start
329  integer, intent(out) :: istrt
330 
331 
332  integer :: len, i
333  real :: lon_strt, tmp(size(lon(:))-1)
334 
335  len = size(lon(:))
336 
337  do i=1,len
338  lon(i) = lon_in_range(lon(i),lon_start)
339  enddo
340 
341  istrt=0
342  do i=1,len-1
343  if (lon(i+1) < lon(i)) then
344  istrt=i+1
345  exit
346  endif
347  enddo
348 
349  if (istrt>1) then ! grid is not monotonic
350  if (abs(lon(len)-lon(1)) < epsln) then
351  tmp = cshift(lon(1:len-1),istrt-1)
352  lon(1:len-1) = tmp
353  lon(len) = lon(1)
354  else
355  lon = cshift(lon,istrt-1)
356  endif
357  lon_strt = lon(1)
358  do i=2,len+1
359  lon(i) = lon_in_range(lon(i),lon_strt)
360  lon_strt = lon(i)
361  enddo
362  endif
363 
364  return
365  end subroutine tranlon
366 
367  function frac_index (value, array)
368  !=======================================================================
369  !
370  ! nearest_index = index of nearest data point within "array" corresponding to
371  ! "value".
372  !
373  ! inputs:
374  !
375  ! value = arbitrary data...same units as elements in "array"
376  ! array = array of data points (must be monotonically increasing)
377  !
378  ! output:
379  !
380  ! nearest_index = index of nearest data point to "value"
381  ! if "value" is outside the domain of "array" then nearest_index = 1
382  ! or "ia" depending on whether array(1) or array(ia) is
383  ! closest to "value"
384  !
385  ! note: if "array" is dimensioned array(0:ia) in the calling
386  ! program, then the returned index should be reduced
387  ! by one to account for the zero base.
388  !
389  ! example:
390  !
391  ! let model depths be defined by the following:
392  ! parameter (km=5)
393  ! dimension z(km)
394  ! data z /5.0, 10.0, 50.0, 100.0, 250.0/
395  !
396  ! k1 = nearest_index (12.5, z, km)
397  ! k2 = nearest_index (0.0, z, km)
398  !
399  ! k1 would be set to 2, and k2 would be set to 1 so that
400  ! z(k1) would be the nearest data point to 12.5 and z(k2) would
401  ! be the nearest data point to 0.0
402  !
403  !=======================================================================
404 
405  integer :: ia, i, ii, unit
406  real :: value, frac_index
407  real, dimension(:) :: array
408  logical keep_going
409 
410  ia = size(array(:))
411 
412  do i=2,ia
413  if (array(i) < array(i-1)) then
414  unit = stdout()
415  write (unit,*) '=> Error: "frac_index" array must be monotonically increasing when searching for nearest value to ',&
416  value
417  write (unit,*) ' array(i) < array(i-1) for i=',i
418  write (unit,*) ' array(i) for i=1..ia follows:'
419  do ii=1,ia
420  write (unit,*) 'i=',ii, ' array(i)=',array(ii)
421  enddo
422  call mpp_error(fatal,' "frac_index" array must be monotonically increasing.')
423  endif
424  enddo
425  if (value < array(1) .or. value > array(ia)) then
426 ! if (value < array(1)) frac_index = 1.
427 ! if (value > array(ia)) frac_index = float(ia)
428  frac_index = -1.0
429  else
430  i=1
431  keep_going = .true.
432  do while (i <= ia .and. keep_going)
433  i = i+1
434  if (value <= array(i)) then
435  frac_index = float(i-1) + (value-array(i-1))/(array(i)-array(i-1))
436  keep_going = .false.
437  endif
438  enddo
439  endif
440  end function frac_index
441 
442  function nearest_index (value, array)
443  !=======================================================================
444  !
445  ! nearest_index = index of nearest data point within "array" corresponding to
446  ! "value".
447  !
448  ! inputs:
449  !
450  ! value = arbitrary data...same units as elements in "array"
451  ! array = array of data points (must be monotonically increasing)
452  ! ia = dimension of "array"
453  !
454  ! output:
455  !
456  ! nearest_index = index of nearest data point to "value"
457  ! if "value" is outside the domain of "array" then nearest_index = 1
458  ! or "ia" depending on whether array(1) or array(ia) is
459  ! closest to "value"
460  !
461  ! note: if "array" is dimensioned array(0:ia) in the calling
462  ! program, then the returned index should be reduced
463  ! by one to account for the zero base.
464  !
465  ! example:
466  !
467  ! let model depths be defined by the following:
468  ! parameter (km=5)
469  ! dimension z(km)
470  ! data z /5.0, 10.0, 50.0, 100.0, 250.0/
471  !
472  ! k1 = nearest_index (12.5, z, km)
473  ! k2 = nearest_index (0.0, z, km)
474  !
475  ! k1 would be set to 2, and k2 would be set to 1 so that
476  ! z(k1) would be the nearest data point to 12.5 and z(k2) would
477  ! be the nearest data point to 0.0
478  !
479  !=======================================================================
480 
481  integer :: nearest_index, ia, i, ii, unit
482  real :: value
483  real, dimension(:) :: array
484  logical keep_going
485 
486  ia = size(array(:))
487 
488  do i=2,ia
489  if (array(i) < array(i-1)) then
490  unit = stdout()
491  write (unit,*) '=> Error: "nearest_index" array must be monotonically increasing &
492  &when searching for nearest value to ',value
493  write (unit,*) ' array(i) < array(i-1) for i=',i
494  write (unit,*) ' array(i) for i=1..ia follows:'
495  do ii=1,ia
496  write (unit,*) 'i=',ii, ' array(i)=',array(ii)
497  enddo
498  call mpp_error(fatal,' "nearest_index" array must be monotonically increasing.')
499  endif
500  enddo
501  if (value < array(1) .or. value > array(ia)) then
502  if (value < array(1)) nearest_index = 1
503  if (value > array(ia)) nearest_index = ia
504  else
505  i=1
506  keep_going = .true.
507  do while (i <= ia .and. keep_going)
508  i = i+1
509  if (value <= array(i)) then
510  nearest_index = i
511  if (array(i)-value > value-array(i-1)) nearest_index = i-1
512  keep_going = .false.
513  endif
514  enddo
515  endif
516  end function nearest_index
517 
518  !#############################################################################
519 
520  subroutine interp_1d_linear(grid1,grid2,data1,data2)
522  real, dimension(:), intent(in) :: grid1, data1, grid2
523  real, dimension(:), intent(inout) :: data2
524 
525  integer :: n1, n2, i, n, ext
526  real :: w
527 
528  n1 = size(grid1(:))
529  n2 = size(grid2(:))
530 
531 
532  do i=2,n1
533  if (grid1(i) <= grid1(i-1)) call mpp_error(fatal, 'grid1 not monotonic')
534  enddo
535 
536  do i=2,n2
537  if (grid2(i) <= grid2(i-1)) call mpp_error(fatal, 'grid2 not monotonic')
538  enddo
539 
540  if (grid1(1) > grid2(1) ) call mpp_error(fatal, 'grid2 lies outside grid1')
541  if (grid1(n1) < grid2(n2) ) call mpp_error(fatal, 'grid2 lies outside grid1')
542 
543  do i=1,n2
544  n = nearest_index(grid2(i),grid1)
545 
546  if (grid1(n) < grid2(i)) then
547  w = (grid2(i)-grid1(n))/(grid1(n+1)-grid1(n))
548  data2(i) = (1.-w)*data1(n) + w*data1(n+1)
549  else
550  if(n==1) then
551  data2(i) = data1(n)
552  else
553  w = (grid2(i)-grid1(n-1))/(grid1(n)-grid1(n-1))
554  data2(i) = (1.-w)*data1(n-1) + w*data1(n)
555  endif
556  endif
557  enddo
558 
559 
560  return
561 
562  end subroutine interp_1d_linear
563 
564  !###################################################################
565  subroutine interp_1d_cubic_spline(grid1, grid2, data1, data2, yp1, ypn)
567  real, dimension(:), intent(in) :: grid1, grid2, data1
568  real, dimension(:), intent(inout) :: data2
569  real, intent(in) :: yp1, ypn
570 
571  real, dimension(size(grid1)) :: y2, u
572  real :: sig, p, qn, un, h, a ,b
573  integer :: n, m, i, k, klo, khi
574 
575  n = size(grid1(:))
576  m = size(grid2(:))
577 
578  do i=2,n
579  if (grid1(i) <= grid1(i-1)) call mpp_error(fatal, 'grid1 not monotonic')
580  enddo
581 
582  do i=2,m
583  if (grid2(i) <= grid2(i-1)) call mpp_error(fatal, 'grid2 not monotonic')
584  enddo
585 
586  if (grid1(1) > grid2(1) ) call mpp_error(fatal, 'grid2 lies outside grid1')
587  if (grid1(n) < grid2(m) ) call mpp_error(fatal, 'grid2 lies outside grid1')
588 
589  if (yp1 >.99e30) then
590  y2(1)=0.
591  u(1)=0.
592  else
593  y2(1)=-0.5
594  u(1)=(3./(grid1(2)-grid1(1)))*((data1(2)-data1(1))/(grid1(2)-grid1(1))-yp1)
595  endif
596 
597  do i=2,n-1
598  sig=(grid1(i)-grid1(i-1))/(grid1(i+1)-grid1(i-1))
599  p=sig*y2(i-1)+2.
600  y2(i)=(sig-1.)/p
601  u(i)=(6.*((data1(i+1)-data1(i))/(grid1(i+1)-grid1(i))-(data1(i)-data1(i-1)) &
602  /(grid1(i)-grid1(i-1)))/(grid1(i+1)-grid1(i-1))-sig*u(i-1))/p
603  enddo
604 
605  if (ypn > .99e30) then
606  qn=0.
607  un=0.
608  else
609  qn=0.5
610  un=(3./(grid1(n)-grid1(n-1)))*(ypn-(data1(n)-data1(n-1))/(grid1(n)-grid1(n-1)))
611  endif
612 
613  y2(n)=(un-qn*u(n-1))/(qn*y2(n-1)+1.)
614 
615  do k=n-1,1,-1
616  y2(k)=y2(k)*y2(k+1)+u(k)
617  enddo
618 
619  do k = 1, m
620  n = nearest_index(grid2(k),grid1)
621  if (grid1(n) < grid2(k)) then
622  klo = n
623  else
624  if(n==1) then
625  klo = n
626  else
627  klo = n -1
628  endif
629  endif
630  khi = klo+1
631  h = grid1(khi)-grid1(klo)
632  a = (grid1(khi) - grid2(k))/h
633  b = (grid2(k) - grid1(klo))/h
634  data2(k) = a*data1(klo) + b*data1(khi)+ ((a**3-a)*y2(klo) + (b**3-b)*y2(khi))*(h**2)/6.
635  enddo
636 
637  end subroutine interp_1d_cubic_spline
638 
639  !###################################################################
640 
641  subroutine interp_1d_1d(grid1,grid2,data1,data2, method, yp1, yp2)
643  real, dimension(:), intent(in) :: grid1, data1, grid2
644  real, dimension(:), intent(inout) :: data2
645  character(len=*), optional, intent(in) :: method
646  real, optional, intent(in) :: yp1, yp2
647 
648  real :: y1, y2
649  character(len=32) :: interp_method
650  integer :: k2, ks, ke
651 
652  k2 = size(grid2(:))
653 
654  interp_method = "linear"
655  if(present(method)) interp_method = method
656  y1 = 1.0e30
657  if(present(yp1)) y1 = yp1
658  y2 = 1.0e30
659  if(present(yp2)) y2 = yp2
660  call find_index(grid1, grid2(1), grid2(k2), ks, ke)
661  select case(trim(interp_method))
662  case("linear")
663  call interp_1d_linear(grid1(ks:ke),grid2,data1(ks:ke),data2)
664  case("cubic_spline")
665  call interp_1d_cubic_spline(grid1(ks:ke),grid2,data1(ks:ke),data2, y1, y2)
666  case default
667  call mpp_error(fatal,"axis_utils: interp_method should be linear or cubic_spline")
668  end select
669 
670  return
671 
672  end subroutine interp_1d_1d
673 
674  !###################################################################
675 
676 
677  subroutine interp_1d_2d(grid1,grid2,data1,data2)
679  real, dimension(:,:), intent(in) :: grid1, data1, grid2
680  real, dimension(:,:), intent(inout) :: data2
681 
682  integer :: n1, n2, i, n, k2, ks, ke
683  real :: w
684 
685  n1 = size(grid1,1)
686  n2 = size(grid2,1)
687  k2 = size(grid2,2)
688 
689  if (n1 /= n2) call mpp_error(fatal,'grid size mismatch')
690 
691  do n=1,n1
692  call find_index(grid1(n,:), grid2(n,1), grid2(n,k2), ks, ke)
693  call interp_1d_linear(grid1(n,ks:ke),grid2(n,:),data1(n,ks:ke),data2(n,:))
694  enddo
695 
696  return
697 
698  end subroutine interp_1d_2d
699 
700  !###################################################################
701 
702  subroutine interp_1d_3d(grid1,grid2,data1,data2, method, yp1, yp2)
704  real, dimension(:,:,:), intent(in) :: grid1, data1, grid2
705  real, dimension(:,:,:), intent(inout) :: data2
706  character(len=*), optional, intent(in) :: method
707  real, optional, intent(in) :: yp1, yp2
708 
709  integer :: n1, n2, m1, m2, k2, i, n, m
710  real :: w, y1, y2
711  character(len=32) :: interp_method
712  integer :: ks, ke
713  n1 = size(grid1,1)
714  n2 = size(grid2,1)
715  m1 = size(grid1,2)
716  m2 = size(grid2,2)
717  k2 = size(grid2,3)
718 
719  interp_method = "linear"
720  if(present(method)) interp_method = method
721  y1 = 1.0e30
722  if(present(yp1)) y1 = yp1
723  y2 = 1.0e30
724  if(present(yp2)) y2 = yp2
725 
726  if (n1 /= n2 .or. m1 /= m2) call mpp_error(fatal,'grid size mismatch')
727 
728  select case(trim(interp_method))
729  case("linear")
730  do m=1,m1
731  do n=1,n1
732  call find_index(grid1(n,m,:), grid2(n,m,1), grid2(n,m,k2), ks, ke)
733  call interp_1d_linear(grid1(n,m,ks:ke),grid2(n,m,:),data1(n,m,ks:ke),data2(n,m,:))
734  enddo
735  enddo
736  case("cubic_spline")
737  do m=1,m1
738  do n=1,n1
739  call find_index(grid1(n,m,:), grid2(n,m,1), grid2(n,m,k2), ks, ke)
740  call interp_1d_cubic_spline(grid1(n,m,ks:ke),grid2(n,m,:), data1(n,m,ks:ke),data2(n,m,:), y1, y2)
741  enddo
742  enddo
743  case default
744  call mpp_error(fatal,"axis_utils: interp_method should be linear or cubic_spline")
745  end select
746 
747  return
748 
749  end subroutine interp_1d_3d
750 
751 
752  !#####################################################################
753  subroutine find_index(grid1, xs, xe, ks, ke)
754  real, dimension(:), intent(in) :: grid1
755  real, intent(in) :: xs, xe
756  integer, intent(out) :: ks, ke
757 
758  integer :: k, nk
759 
760  nk = size(grid1(:))
761 
762  ks = 0; ke = 0
763  do k = 1, nk-1
764  if(grid1(k) <= xs .and. grid1(k+1) > xs ) then
765  ks = k
766  exit
767  endif
768  enddo
769  do k = nk, 2, -1
770  if(grid1(k) >= xe .and. grid1(k-1) < xe ) then
771  ke = k
772  exit
773  endif
774  enddo
775 
776  if(ks == 0 ) call mpp_error(fatal,' xs locate outside of grid1')
777  if(ke == 0 ) call mpp_error(fatal,' xe locate outside of grid1')
778 
779  end subroutine find_index
780 
781 end module axis_utils_mod
782 
783 #ifdef test_axis_utils
784 
785 program test
786 
787 use fms_mod, only : fms_init, file_exist, open_namelist_file, check_nml_error
788 use fms_mod, only : close_file
789 use mpp_mod, only : mpp_error, fatal, stdout
790 use mpp_mod, only : input_nml_file
791 use axis_utils_mod, only: interp_1d
792 
793 implicit none
794 
795 
796 
797 integer, parameter :: maxsize = 100
798 
799 integer :: n_src = 0
800 integer :: n_dst = 0
801 real, dimension(MAXSIZE) :: grid_src = 0
802 real, dimension(MAXSIZE) :: grid_dst = 0
803 real, dimension(MAXSIZE) :: data_src = 0
804 
805 namelist / test_axis_utils_nml / n_src, n_dst, grid_src, grid_dst, data_src
806 
807 real, allocatable :: data_dst(:)
808 integer :: unit, ierr, io
809 
810  call fms_init()
811 
812  !--- default option of data
813  n_src = 31
814  n_dst = 40
815  grid_src(1:n_src) = (/ -63.6711465476916, -63.6711455476916, 166.564180735096, 401.25299580552, &
816  641.056493022762, 886.219516665347, 1137.35352761133, 1394.4936854079, &
817  1657.17893448689, 1925.64572676068, 2200.13183483549, 2480.9124139255, &
818  2768.35396680912, 3062.86513953019, 3675.47369643284, 4325.10564183322, &
819  5020.19039479527, 5769.85432323481, 6584.25101514851, 7475.94655633703, &
820  8462.01951335773, 9568.28246037887, 10178.3869413515, 10834.1425668942, &
821  11543.5265942777, 12317.3907407535, 13170.4562394288, 14125.6466646843, &
822  15225.8720618086, 16554.7859690842, 19697.1334102613 /)
823  grid_dst(1:n_dst) = (/ 1002.9522552602, 1077.51144617887, 1163.37842788755, 1264.19848463606, &
824  1382.57557953916, 1521.56713587855, 1684.76300370633, 1876.37817787584, &
825  2101.36166220498, 2365.52429149707, 2675.68881278444, 3039.86610206727, &
826  3467.4620678435, 3969.52058529847, 4553.81573511231, 5159.54844211827, &
827  5765.28114912423, 6371.01385613019, 6976.74656313614, 7582.4792701421, &
828  8188.21197714806, 8793.94468415402, 9399.67739115997, 10005.4100981659, &
829  10611.1428051719, 11216.8755121778, 11822.6082191838, 12428.3409261898, &
830  13034.0736331957, 13639.8063402017, 14245.5390472076, 14851.2717542136, &
831  15457.0044612196, 16062.7371682255, 16668.4698752315, 17274.2025822374, &
832  17879.9352892434, 18485.6679962493, 19091.4007032553, 19697.1334102613 /)
833  data_src(1:n_src) = (/ 309.895999643929, 309.991081541887, 309.971074746584, 310.873654697145, &
834  311.946530606618, 312.862249229647, 314.821236806913, 315.001269608758, &
835  315.092410930288, 315.19010999336, 315.122964496815, 315.057882573487, &
836  314.998796850493, 314.984586411292, 315.782246062002, 318.142544345795, &
837  321.553905292867, 325.247730854554, 329.151282227113, 332.835673638378, &
838  336.810414210932, 341.64530983048, 344.155248759994, 346.650476976385, &
839  349.106430095269, 351.915323032738, 354.709396583792, 359.68904432446, &
840  371.054289820675, 395.098187506342, 446.150726850039 /)
841 
842 
843  !---reading namelist
844 #ifdef INTERNAL_FILE_NML
845  read (input_nml_file, test_axis_utils_nml, iostat=io)
846  ierr = check_nml_error(io,'test_axis_utils_nml')
847 #else
848  if(file_exist('input.nml')) then
849  unit = open_namelist_file()
850  ierr=1
851  do while (ierr /= 0)
852  read (unit, nml=test_axis_utils_nml, iostat=io, end=10)
853  ierr = check_nml_error(io,'test_axis_utils_nml') ! also initializes nml error codes
854  enddo
855  10 call close_file(unit)
856  endif
857 #endif
858 
859  if(n_src >maxsize) call mpp_error(fatal, 'test_axis_utils: nml n_src is greater than MAXSIZE')
860  if(n_dst >maxsize) call mpp_error(fatal, 'test_axis_utils: nml n_dst is greater than MAXSIZE')
861 
862  allocate(data_dst(n_dst) )
863 
864 
865 
866  !--- write out data
867  unit = stdout()
868  write(unit,*)' the source grid is ', grid_src(1:n_src)
869  write(unit,*)' the destination grid is ', grid_dst(1:n_dst)
870  write(unit,*)' the source data is ', data_src(1:n_src)
871  call interp_1d(grid_src(1:n_src), grid_dst(1:n_dst), data_src(1:n_src), data_dst, "linear")
872  write(unit,*)' the destination data using linear interpolation is ', data_dst(1:n_dst)
873  call interp_1d(grid_src(1:n_src), grid_dst(1:n_dst), data_src(1:n_src), data_dst, "cubic_spline")
874  write(unit,*)' the destination data using cublic spline interpolation is ', data_dst(1:n_dst)
875 
876 end program test
877 
878 
879 #endif /* test_axis_utils */
880 
881 
882 
883 
Definition: fms.F90:20
type(atttype), save, public default_att
Definition: mpp_io.F90:1073
integer, parameter maxatts
Definition: axis_utils.F90:62
integer function, public nearest_index(value, array)
Definition: axis_utils.F90:443
logical function, public get_axis_modulo(axis)
Definition: axis_utils.F90:203
logical function, public fms_error_handler(routine, message, err_msg)
Definition: fms.F90:573
logical function, public get_axis_modulo_times(axis, tbeg, tend)
Definition: axis_utils.F90:225
subroutine interp_1d_1d(grid1, grid2, data1, data2, method, yp1, yp2)
Definition: axis_utils.F90:642
Definition: mpp.F90:39
integer function, public check_nml_error(IOSTAT, NML_NAME)
Definition: fms.F90:658
subroutine interp_1d_3d(grid1, grid2, data1, data2, method, yp1, yp2)
Definition: axis_utils.F90:703
character(len=input_str_length), dimension(:), allocatable, target, public input_nml_file
Definition: mpp.F90:1378
real, parameter epsln
Definition: axis_utils.F90:63
subroutine, public fms_init(localcomm)
Definition: fms.F90:353
subroutine find_index(grid1, xs, xe, ks, ke)
Definition: axis_utils.F90:754
subroutine interp_1d_cubic_spline(grid1, grid2, data1, data2, yp1, ypn)
Definition: axis_utils.F90:566
subroutine, public get_axis_cart(axis, cart)
Definition: axis_utils.F90:79
type(axistype), save, public default_axis
Definition: mpp_io.F90:1071
logical function, public get_axis_fold(axis)
Definition: axis_utils.F90:269
subroutine interp_1d_linear(grid1, grid2, data1, data2)
Definition: axis_utils.F90:521
subroutine interp_1d_2d(grid1, grid2, data1, data2)
Definition: axis_utils.F90:678
************************************************************************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)
subroutine, public get_axis_bounds(axis, axis_bound, axes, bnd_name, err_msg)
Definition: axis_utils.F90:149
real function, public frac_index(value, array)
Definition: axis_utils.F90:368
real, parameter fp5
Definition: axis_utils.F90:64
logical function, public string_array_index(string, string_array, index)
Definition: fms.F90:831
real, parameter f360
Definition: axis_utils.F90:64
real function, public lon_in_range(lon, l_strt)
Definition: axis_utils.F90:291
subroutine, public tranlon(lon, lon_start, istrt)
Definition: axis_utils.F90:319