30#ifdef use_deprecated_io
36 use mpp_mod,
only:
mpp_error, fatal, stdout
37 use fms_mod,
only: lowercase, string_array_index, fms_error_handler
41 public get_axis_cart, get_axis_bounds, get_axis_modulo, get_axis_fold, lon_in_range, &
42 tranlon, frac_index, nearest_index, interp_1d, get_axis_modulo_times
46 integer,
parameter :: maxatts = 100
47 real,
parameter :: epsln= 1.e-10
48 real,
parameter :: fp5 = 0.5, f360 = 360.0
51#include<file_version.h>
65 module procedure interp_1d_1d
66 module procedure interp_1d_2d
67 module procedure interp_1d_3d
75 subroutine get_axis_cart(axis, cart)
77 type(axistype),
intent(in) :: axis
78 character(len=1),
intent(out) :: cart
79 character(len=1) :: axis_cart
80 character(len=16),
dimension(2) :: lon_names, lat_names
81 character(len=16),
dimension(3) :: z_names
82 character(len=16),
dimension(2) :: t_names
83 character(len=16),
dimension(3) :: lon_units, lat_units
84 character(len=8) ,
dimension(4) :: z_units
85 character(len=3) ,
dimension(6) :: t_units
86 character(len=32) :: name
89 lon_names = (/
'lon',
'x '/)
90 lat_names = (/
'lat',
'y '/)
91 z_names = (/
'depth ',
'height',
'z '/)
92 t_names = (/
'time',
't '/)
93 lon_units = (/
'degrees_e ',
'degrees_east',
'degreese '/)
94 lat_units = (/
'degrees_n ',
'degrees_north',
'degreesn '/)
95 z_units = (/
'cm ',
'm ',
'pa ',
'hpa'/)
96 t_units = (/
'sec',
'min',
'hou',
'day',
'mon',
'yea'/)
101 if ( lowercase(axis_cart) ==
'x' ) cart =
'X'
102 if ( lowercase(axis_cart) ==
'y' ) cart =
'Y'
103 if ( lowercase(axis_cart) ==
'z' ) cart =
'Z'
104 if ( lowercase(axis_cart) ==
't' ) cart =
'T'
106 if (cart /=
'X' .and. cart /=
'Y' .and. cart /=
'Z' .and. cart /=
'T')
then
108 name = lowercase(name)
109 do i=1,
size(lon_names(:))
110 if (trim(name(1:3)) == trim(lon_names(i))) cart =
'X'
112 do i=1,
size(lat_names(:))
113 if (trim(name(1:3)) == trim(lat_names(i))) cart =
'Y'
115 do i=1,
size(z_names(:))
116 if (trim(name) == trim(z_names(i))) cart =
'Z'
118 do i=1,
size(t_names(:))
119 if (trim(name) == t_names(i)) cart =
'T'
123 if (cart /=
'X' .and. cart /=
'Y' .and. cart /=
'Z' .and. cart /=
'T')
then
125 name = lowercase(name)
126 do i=1,
size(lon_units(:))
127 if (trim(name) == trim(lon_units(i))) cart =
'X'
129 do i=1,
size(lat_units(:))
130 if (trim(name) == trim(lat_units(i))) cart =
'Y'
132 do i=1,
size(z_units(:))
133 if (trim(name) == trim(z_units(i))) cart =
'Z'
135 do i=1,
size(t_units(:))
136 if (name(1:3) == trim(t_units(i))) cart =
'T'
142 end subroutine get_axis_cart
145 subroutine get_axis_bounds(axis,axis_bound,axes,bnd_name,err_msg)
147 type(axistype),
intent(in) :: axis
148 type(axistype),
intent(inout) :: axis_bound
149 type(axistype),
intent(in),
dimension(:) :: axes
150 character(len=*),
intent(inout),
optional :: bnd_name
151 character(len=*),
intent(out),
optional :: err_msg
153 real,
dimension(:),
allocatable :: data, tmp
156 character(len=128) :: name, units
157 character(len=256) :: longname
158 character(len=1) :: cartesian
159 logical :: bounds_found
161 if(
present(err_msg))
then
164 axis_bound = default_axis
166 cartesian=cartesian, len=len)
167 if(len .LE. 0)
return
168 allocate(
data(len+1))
170 bounds_found = mpp_get_axis_bounds(axis,
data, name=name)
171 longname = trim(longname)//
' bounds'
173 if(.not.bounds_found .and. len>1 )
then
176 name = trim(name)//
'_bnds'
178 call mpp_get_axis_data(axis,tmp)
180 data(i)= tmp(i-1)+fp5*(tmp(i)-tmp(i-1))
182 data(1)= tmp(1)- fp5*(tmp(2)-tmp(1))
183 if (abs(
data(1)) < epsln)
data(1) = 0.0
184 data(len+1)= tmp(len)+ fp5*(tmp(len)-tmp(len-1))
185 if (
data(1) == 0.0)
then
186 if (abs(
data(len+1)-360.) > epsln)
data(len+1)=360.0
189 if(bounds_found .OR. len>1)
then
191 longname,cartesian=cartesian,data=data)
193 if(
allocated(tmp))
deallocate(tmp)
197 end subroutine get_axis_bounds
202 function get_axis_modulo(axis)
204 type(axistype) :: axis
205 logical :: get_axis_modulo
207 type(atttype),
dimension(:),
allocatable :: atts
214 get_axis_modulo=.false.
216 if (lowercase(trim(
mpp_get_att_name(atts(i)))) ==
'modulo') get_axis_modulo = .true.
222 end function get_axis_modulo
225 function get_axis_modulo_times(axis, tbeg, tend)
227 logical :: get_axis_modulo_times
228 type(axistype),
intent(in) :: axis
229 character(len=*),
intent(out) :: tbeg, tend
231 type(atttype),
dimension(:),
allocatable :: atts
232 logical :: found_tbeg, found_tend
244 call mpp_error(fatal,
'error in get: len(tbeg) too small to hold attribute')
251 call mpp_error(fatal,
'error in get: len(tend) too small to hold attribute')
258 if(found_tbeg .and. .not.found_tend)
then
259 call mpp_error(fatal,
'error in get: Found modulo_beg but not modulo_end')
261 if(.not.found_tbeg .and. found_tend)
then
262 call mpp_error(fatal,
'error in get: Found modulo_end but not modulo_beg')
265 get_axis_modulo_times = found_tbeg
267 end function get_axis_modulo_times
271 function get_axis_fold(axis)
273 type(axistype) :: axis
274 logical :: get_axis_fold
276 type(atttype),
dimension(:),
allocatable :: atts
283 get_axis_fold=.false.
291 end function get_axis_fold
295 function lon_in_range(lon, l_strt)
296 real :: lon, l_strt, lon_in_range, l_end
301 if (abs(lon_in_range - l_strt) < 1.e-4)
then
302 lon_in_range = l_strt
306 if (abs(lon_in_range - l_end) < 1.e-4)
then
307 lon_in_range = l_strt
312 if (lon_in_range < l_strt)
then
313 lon_in_range = lon_in_range + f360;
314 else if (lon_in_range > l_end)
then
315 lon_in_range = lon_in_range - f360;
321 end function lon_in_range
324 subroutine tranlon(lon, lon_start, istrt)
333 real,
intent(inout),
dimension(:) :: lon
334 real,
intent(in) :: lon_start
335 integer,
intent(out) :: istrt
339 real :: lon_strt, tmp(size(lon(:))-1)
344 lon(i) = lon_in_range(lon(i),lon_start)
349 if (lon(i+1) < lon(i))
then
356 if (abs(lon(len)-lon(1)) < epsln)
then
357 tmp = cshift(lon(1:len-1),istrt-1)
361 lon = cshift(lon,istrt-1)
365 lon(i) = lon_in_range(lon(i),lon_strt)
371 end subroutine tranlon
407 function frac_index (value, array)
408 integer :: ia, i, ii, unit
411 real,
dimension(:) :: array
417 if (array(i) < array(i-1))
then
420 '=> Error: "frac_index" array must be monotonically increasing when searching for nearest value to ',
value
421 write (unit,*)
' array(i) < array(i-1) for i=',i
422 write (unit,*)
' array(i) for i=1..ia follows:'
424 write (unit,*)
'i=',ii,
' array(i)=',array(ii)
426 call mpp_error(fatal,
' "frac_index" array must be monotonically increasing.')
429 if (
value < array(1) .or.
value > array(ia))
then
436 do while (i <= ia .and. keep_going)
438 if (
value <= array(i))
then
439 frac_index = float(i-1) + (
value-array(i-1))/(array(i)-array(i-1))
444 end function frac_index
482 function nearest_index (value, array)
488 integer :: nearest_index
492 Integer :: i, ii, unit
495 real,
dimension(:) :: array
501 if (array(i) < array(i-1))
then
503 write (unit,*)
'=> Error: "nearest_index" array must be monotonically increasing &
504 &when searching for nearest value to ',
value
505 write (unit,*)
' array(i) < array(i-1) for i=',i
506 write (unit,*)
' array(i) for i=1..ia follows:'
508 write (unit,*)
'i=',ii,
' array(i)=',array(ii)
510 call mpp_error(fatal,
' "nearest_index" array must be monotonically increasing.')
513 if (
value < array(1) .or.
value > array(ia))
then
514 if (
value < array(1)) nearest_index = 1
515 if (
value > array(ia)) nearest_index = ia
519 do while (i <= ia .and. keep_going)
521 if (
value <= array(i))
then
523 if (array(i)-
value >
value-array(i-1)) nearest_index = i-1
528 end function nearest_index
532 subroutine interp_1d_linear(grid1,grid2,data1,data2)
534 real,
dimension(:),
intent(in) :: grid1, data1, grid2
535 real,
dimension(:),
intent(inout) :: data2
537 integer :: n1, n2, i, n
545 if (grid1(i) <= grid1(i-1))
call mpp_error(fatal,
'grid1 not monotonic')
549 if (grid2(i) <= grid2(i-1))
call mpp_error(fatal,
'grid2 not monotonic')
552 if (grid1(1) > grid2(1) )
call mpp_error(fatal,
'grid2 lies outside grid1')
553 if (grid1(n1) < grid2(n2) )
call mpp_error(fatal,
'grid2 lies outside grid1')
556 n = nearest_index(grid2(i),grid1)
558 if (grid1(n) < grid2(i))
then
559 w = (grid2(i)-grid1(n))/(grid1(n+1)-grid1(n))
560 data2(i) = (1.-w)*data1(n) + w*data1(n+1)
565 w = (grid2(i)-grid1(n-1))/(grid1(n)-grid1(n-1))
566 data2(i) = (1.-w)*data1(n-1) + w*data1(n)
574 end subroutine interp_1d_linear
577 subroutine interp_1d_cubic_spline(grid1, grid2, data1, data2, yp1, ypn)
579 real,
dimension(:),
intent(in) :: grid1, grid2, data1
580 real,
dimension(:),
intent(inout) :: data2
581 real,
intent(in) :: yp1, ypn
583 real,
dimension(size(grid1)) :: y2, u
584 real :: sig, p, qn, un, h, a ,b
585 integer :: n, m, i, k, klo, khi
591 if (grid1(i) <= grid1(i-1))
call mpp_error(fatal,
'grid1 not monotonic')
595 if (grid2(i) <= grid2(i-1))
call mpp_error(fatal,
'grid2 not monotonic')
598 if (grid1(1) > grid2(1) )
call mpp_error(fatal,
'grid2 lies outside grid1')
599 if (grid1(n) < grid2(m) )
call mpp_error(fatal,
'grid2 lies outside grid1')
601 if (yp1 >.99e30)
then
606 u(1)=(3./(grid1(2)-grid1(1)))*((data1(2)-data1(1))/(grid1(2)-grid1(1))-yp1)
610 sig=(grid1(i)-grid1(i-1))/(grid1(i+1)-grid1(i-1))
613 u(i)=(6.*((data1(i+1)-data1(i))/(grid1(i+1)-grid1(i))-(data1(i)-data1(i-1)) &
614 /(grid1(i)-grid1(i-1)))/(grid1(i+1)-grid1(i-1))-sig*u(i-1))/p
617 if (ypn > .99e30)
then
622 un=(3./(grid1(n)-grid1(n-1)))*(ypn-(data1(n)-data1(n-1))/(grid1(n)-grid1(n-1)))
625 y2(n)=(un-qn*u(n-1))/(qn*y2(n-1)+1.)
628 y2(k)=y2(k)*y2(k+1)+u(k)
632 n = nearest_index(grid2(k),grid1)
633 if (grid1(n) < grid2(k))
then
643 h = grid1(khi)-grid1(klo)
644 a = (grid1(khi) - grid2(k))/h
645 b = (grid2(k) - grid1(klo))/h
646 data2(k) = a*data1(klo) + b*data1(khi)+ ((a**3-a)*y2(klo) + (b**3-b)*y2(khi))*(h**2)/6.
649 end subroutine interp_1d_cubic_spline
653 subroutine interp_1d_1d(grid1,grid2,data1,data2, method, yp1, yp2)
655 real,
dimension(:),
intent(in) :: grid1, data1, grid2
656 real,
dimension(:),
intent(inout) :: data2
657 character(len=*),
optional,
intent(in) :: method
658 real,
optional,
intent(in) :: yp1, yp2
661 character(len=32) :: interp_method
662 integer :: k2, ks, ke
666 interp_method =
"linear"
667 if(
present(method)) interp_method = method
669 if(
present(yp1)) y1 = yp1
671 if(
present(yp2)) y2 = yp2
672 call find_index(grid1, grid2(1), grid2(k2), ks, ke)
673 select case(trim(interp_method))
675 call interp_1d_linear(grid1(ks:ke),grid2,data1(ks:ke),data2)
677 call interp_1d_cubic_spline(grid1(ks:ke),grid2,data1(ks:ke),data2, y1, y2)
679 call mpp_error(fatal,
"axis_utils: interp_method should be linear or cubic_spline")
684 end subroutine interp_1d_1d
689 subroutine interp_1d_2d(grid1,grid2,data1,data2)
691 real,
dimension(:,:),
intent(in) :: grid1, data1, grid2
692 real,
dimension(:,:),
intent(inout) :: data2
694 integer :: n1, n2, n, k2, ks, ke
700 if (n1 /= n2)
call mpp_error(fatal,
'grid size mismatch')
703 call find_index(grid1(n,:), grid2(n,1), grid2(n,k2), ks, ke)
704 call interp_1d_linear(grid1(n,ks:ke),grid2(n,:),data1(n,ks:ke),data2(n,:))
709 end subroutine interp_1d_2d
713 subroutine interp_1d_3d(grid1,grid2,data1,data2, method, yp1, yp2)
715 real,
dimension(:,:,:),
intent(in) :: grid1, data1, grid2
716 real,
dimension(:,:,:),
intent(inout) :: data2
717 character(len=*),
optional,
intent(in) :: method
718 real,
optional,
intent(in) :: yp1, yp2
720 integer :: n1, n2, m1, m2, k2, n, m
722 character(len=32) :: interp_method
730 interp_method =
"linear"
731 if(
present(method)) interp_method = method
733 if(
present(yp1)) y1 = yp1
735 if(
present(yp2)) y2 = yp2
737 if (n1 /= n2 .or. m1 /= m2)
call mpp_error(fatal,
'grid size mismatch')
739 select case(trim(interp_method))
743 call find_index(grid1(n,m,:), grid2(n,m,1), grid2(n,m,k2), ks, ke)
744 call interp_1d_linear(grid1(n,m,ks:ke),grid2(n,m,:),data1(n,m,ks:ke),data2(n,m,:))
750 call find_index(grid1(n,m,:), grid2(n,m,1), grid2(n,m,k2), ks, ke)
751 call interp_1d_cubic_spline(grid1(n,m,ks:ke),grid2(n,m,:), data1(n,m,ks:ke),data2(n,m,:), y1, y2)
755 call mpp_error(fatal,
"axis_utils: interp_method should be linear or cubic_spline")
760 end subroutine interp_1d_3d
764 subroutine find_index(grid1, xs, xe, ks, ke)
765 real,
dimension(:),
intent(in) :: grid1
766 real,
intent(in) :: xs, xe
767 integer,
intent(out) :: ks, ke
775 if(grid1(k) <= xs .and. grid1(k+1) > xs )
then
781 if(grid1(k) >= xe .and. grid1(k-1) < xe )
then
787 if(ks == 0 )
call mpp_error(fatal,
' xs locate outside of grid1')
788 if(ke == 0 )
call mpp_error(fatal,
' xe locate outside of grid1')
790 end subroutine find_index
792end module axis_utils_mod
integer function mpp_get_att_length(att)
return the length of an attribute.
character(len=len(att%name)) function mpp_get_att_name(att)
return the name of an attribute.
integer function mpp_get_att_type(att)
return the type of an attribute.
character(len=att%len) function mpp_get_att_char(att)
return the char value of an attribute.
Get file global metadata.