FMS 2025.01-dev
Flexible Modeling System
Loading...
Searching...
No Matches
axis_utils.F90
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!> @defgroup axis_utils_mod axis_utils_mod
20!> @ingroup axis_utils
21!> @brief A set of utilities for manipulating axes and extracting axis attributes,
22!! A more recent version of this module using the updated fms2_io is available, @ref axis_utils2
23!! but this module should be used if still using @ref mpp_io
24!!
25!> @author M.J. Harrison
26
27!> @addtogroup axis_utils_mod
28!> @{
29module axis_utils_mod
30#ifdef use_deprecated_io
31 use netcdf
32 use mpp_io_mod, only: axistype, atttype, default_axis, default_att, &
33 mpp_get_atts, mpp_get_axis_data, mpp_modify_meta, &
35 mpp_get_att_length, mpp_get_axis_bounds
36 use mpp_mod, only: mpp_error, fatal, stdout
37 use fms_mod, only: lowercase, string_array_index, fms_error_handler
38
39 implicit none
40
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
43
44 private
45
46 integer, parameter :: maxatts = 100
47 real, parameter :: epsln= 1.e-10
48 real, parameter :: fp5 = 0.5, f360 = 360.0
49
50! Include variable "version" to be written to log file.
51#include<file_version.h>
52
53!> @}
54
55 !> Perform 1D interpolation between grids.
56 !!
57 !> Data and grids can have 1, 2, or 3 dimensions.
58 !! @param grid1
59 !! @param grid2
60 !! @param data1 Data to interpolate
61 !! @param [inout] data2 Interpolated data
62 !! @param method Either "linear" or "cubic_spline" interpolation method, default="linear"
63 !! @ingroup axis_utils_mod
64 interface interp_1d
65 module procedure interp_1d_1d
66 module procedure interp_1d_2d
67 module procedure interp_1d_3d
68 end interface
69!> @addtogroup axis_utils_mod
70!> @{
71
72contains
73
74 !> @brief Returns X,Y,Z or T cartesian attribute
75 subroutine get_axis_cart(axis, cart)
76
77 type(axistype), intent(in) :: axis !< axis to get data from
78 character(len=1), intent(out) :: cart !< Returned cartesian axis
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
87 integer :: i
88
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'/)
97
98 call mpp_get_atts(axis,cartesian=axis_cart)
99 cart = 'N'
100
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'
105
106 if (cart /= 'X' .and. cart /= 'Y' .and. cart /= 'Z' .and. cart /= 'T') then
107 call mpp_get_atts(axis,name=name)
108 name = lowercase(name)
109 do i=1,size(lon_names(:))
110 if (trim(name(1:3)) == trim(lon_names(i))) cart = 'X'
111 enddo
112 do i=1,size(lat_names(:))
113 if (trim(name(1:3)) == trim(lat_names(i))) cart = 'Y'
114 enddo
115 do i=1,size(z_names(:))
116 if (trim(name) == trim(z_names(i))) cart = 'Z'
117 enddo
118 do i=1,size(t_names(:))
119 if (trim(name) == t_names(i)) cart = 'T'
120 enddo
121 end if
122
123 if (cart /= 'X' .and. cart /= 'Y' .and. cart /= 'Z' .and. cart /= 'T') then
124 call mpp_get_atts(axis,units=name)
125 name = lowercase(name)
126 do i=1,size(lon_units(:))
127 if (trim(name) == trim(lon_units(i))) cart = 'X'
128 enddo
129 do i=1,size(lat_units(:))
130 if (trim(name) == trim(lat_units(i))) cart = 'Y'
131 enddo
132 do i=1,size(z_units(:))
133 if (trim(name) == trim(z_units(i))) cart = 'Z'
134 enddo
135 do i=1,size(t_units(:))
136 if (name(1:3) == trim(t_units(i))) cart = 'T'
137 enddo
138 end if
139
140 return
141
142 end subroutine get_axis_cart
143
144 !> @brief Return axis_bound either from an array of available axes, or defined based on axis mid-points
145 subroutine get_axis_bounds(axis,axis_bound,axes,bnd_name,err_msg)
146
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
152
153 real, dimension(:), allocatable :: data, tmp
154
155 integer :: i, len
156 character(len=128) :: name, units
157 character(len=256) :: longname
158 character(len=1) :: cartesian
159 logical :: bounds_found
160
161 if(present(err_msg)) then
162 err_msg = ''
163 endif
164 axis_bound = default_axis
165 call mpp_get_atts(axis,units=units,longname=longname,&
166 cartesian=cartesian, len=len)
167 if(len .LE. 0) return
168 allocate(data(len+1))
169
170 bounds_found = mpp_get_axis_bounds(axis, data, name=name)
171 longname = trim(longname)//' bounds'
172
173 if(.not.bounds_found .and. len>1 ) then
174 ! The following calculation can not be done for len=1
175 call mpp_get_atts(axis,name=name)
176 name = trim(name)//'_bnds'
177 allocate(tmp(len))
178 call mpp_get_axis_data(axis,tmp)
179 do i=2,len
180 data(i)= tmp(i-1)+fp5*(tmp(i)-tmp(i-1))
181 enddo
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
187 endif
188 endif
189 if(bounds_found .OR. len>1) then
190 call mpp_modify_meta(axis_bound,name=name,units=units,longname=&
191 longname,cartesian=cartesian,data=data)
192 endif
193 if(allocated(tmp)) deallocate(tmp)
194 deallocate(data)
195
196 return
197 end subroutine get_axis_bounds
198
199 !> @brief Checks if 'modulo' variable exists for a given axis.
200 !!
201 !> @return true if modulo variable exists in fileobj for the given axis name.
202 function get_axis_modulo(axis)
203
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 !> @return logical get_axis_modulo_times
225 function get_axis_modulo_times(axis, tbeg, tend)
226
227 logical :: get_axis_modulo_times
228 type(axistype), intent(in) :: axis
229 character(len=*), intent(out) :: tbeg, tend
230 integer :: natt, i
231 type(atttype), dimension(:), allocatable :: atts
232 logical :: found_tbeg, found_tend
233
234 call mpp_get_atts(axis,natts=natt)
235 allocate(atts(natt))
236 call mpp_get_atts(axis,atts=atts)
237
238 found_tbeg = .false.
239 found_tend = .false.
240
241 do i = 1,natt
242 if(lowercase(trim(mpp_get_att_name(atts(i)))) == 'modulo_beg') then
243 if(mpp_get_att_length(atts(i)) > len(tbeg)) then
244 call mpp_error(fatal,'error in get: len(tbeg) too small to hold attribute')
245 endif
246 tbeg = trim(mpp_get_att_char(atts(i)))
247 found_tbeg = .true.
248 endif
249 if(lowercase(trim(mpp_get_att_name(atts(i)))) == 'modulo_end') then
250 if(mpp_get_att_length(atts(i)) > len(tend)) then
251 call mpp_error(fatal,'error in get: len(tend) too small to hold attribute')
252 endif
253 tend = trim(mpp_get_att_char(atts(i)))
254 found_tend = .true.
255 endif
256 enddo
257
258 if(found_tbeg .and. .not.found_tend) then
259 call mpp_error(fatal,'error in get: Found modulo_beg but not modulo_end')
260 endif
261 if(.not.found_tbeg .and. found_tend) then
262 call mpp_error(fatal,'error in get: Found modulo_end but not modulo_beg')
263 endif
264
265 get_axis_modulo_times = found_tbeg
266
267 end function get_axis_modulo_times
268
269 !> @brief Returns if axis is folded at a boundary (non-standard meta-data)
270 !! @return logical get_axis_fold
271 function get_axis_fold(axis)
272
273 type(axistype) :: axis
274 logical :: get_axis_fold
275 integer :: natt, i
276 type(atttype), dimension(:), allocatable :: atts
277
278
279 call mpp_get_atts(axis,natts=natt)
280 allocate(atts(natt))
281 call mpp_get_atts(axis,atts=atts)
282
283 get_axis_fold=.false.
284 do i = 1,natt
285 if (mpp_get_att_char(atts(i)) == 'fold_top') get_axis_fold = .true.
286 enddo
287
288 deallocate(atts)
289
290 return
291 end function get_axis_fold
292
293 !> @brief Returns lon_strt <= longitude <= lon_strt+360
294 !! @return real lon_in_range
295 function lon_in_range(lon, l_strt)
296 real :: lon, l_strt, lon_in_range, l_end
297
298 lon_in_range = lon
299 l_end = l_strt+360.
300
301 if (abs(lon_in_range - l_strt) < 1.e-4) then
302 lon_in_range = l_strt
303 return
304 endif
305
306 if (abs(lon_in_range - l_end) < 1.e-4) then
307 lon_in_range = l_strt
308 return
309 endif
310
311 do
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;
316 else
317 exit
318 end if
319 end do
320
321 end function lon_in_range
322
323 !> @brief Returns monotonic array of longitudes s.t., lon_strt <= lon(:) <= lon_strt+360.
324 subroutine tranlon(lon, lon_start, istrt)
325
326 ! returns array of longitudes s.t. lon_strt <= lon < lon_strt+360.
327 ! also, the first istrt-1 entries are moved to the end of the array
328 !
329 ! e.g.
330 ! lon = 0 1 2 3 4 5 ... 358 359; lon_strt = 3 ==>
331 ! tranlon = 3 4 5 6 7 8 ... 359 360 361 362; istrt = 4
332
333 real, intent(inout), dimension(:) :: lon
334 real, intent(in) :: lon_start
335 integer, intent(out) :: istrt
336
337
338 integer :: len, i
339 real :: lon_strt, tmp(size(lon(:))-1)
340
341 len = size(lon(:))
342
343 do i=1,len
344 lon(i) = lon_in_range(lon(i),lon_start)
345 enddo
346
347 istrt=0
348 do i=1,len-1
349 if (lon(i+1) < lon(i)) then
350 istrt=i+1
351 exit
352 endif
353 enddo
354
355 if (istrt>1) then ! grid is not monotonic
356 if (abs(lon(len)-lon(1)) < epsln) then
357 tmp = cshift(lon(1:len-1),istrt-1)
358 lon(1:len-1) = tmp
359 lon(len) = lon(1)
360 else
361 lon = cshift(lon,istrt-1)
362 endif
363 lon_strt = lon(1)
364 do i=2,len+1
365 lon(i) = lon_in_range(lon(i),lon_strt)
366 lon_strt = lon(i)
367 enddo
368 endif
369
370 return
371 end subroutine tranlon
372
373 !> nearest_index = index of nearest data point within "array" corresponding to
374 !! "value".
375 !!
376 !! inputs:
377 !!
378 !! value = arbitrary data...same units as elements in "array"
379 !! array = array of data points (must be monotonically increasing)
380 !!
381 !! output:
382 !!
383 !! nearest_index = index of nearest data point to "value"
384 !! if "value" is outside the domain of "array" then nearest_index = 1
385 !! or "ia" depending on whether array(1) or array(ia) is
386 !! closest to "value"
387 !!
388 !! note: if "array" is dimensioned array(0:ia) in the calling
389 !! program, then the returned index should be reduced
390 !! by one to account for the zero base.
391 !!
392 !! example:
393 !!
394 !! let model depths be defined by the following:
395 !! parameter (km=5)
396 !! dimension z(km)
397 !! data z /5.0, 10.0, 50.0, 100.0, 250.0/
398 !!
399 !! k1 = nearest_index (12.5, z, km)
400 !! k2 = nearest_index (0.0, z, km)
401 !!
402 !! k1 would be set to 2, and k2 would be set to 1 so that
403 !! z(k1) would be the nearest data point to 12.5 and z(k2) would
404 !! be the nearest data point to 0.0
405 !!
406 !! @return real frac_index
407 function frac_index (value, array)
408 integer :: ia, i, ii, unit
409 real :: value !< arbitrary data...same units as elements in "array"
410 real :: frac_index
411 real, dimension(:) :: array !< array of data points (must be monotonically increasing)
412 logical keep_going
413
414 ia = size(array(:))
415
416 do i=2,ia
417 if (array(i) < array(i-1)) then
418 unit = stdout()
419 write (unit,*) &
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:'
423 do ii=1,ia
424 write (unit,*) 'i=',ii, ' array(i)=',array(ii)
425 enddo
426 call mpp_error(fatal,' "frac_index" array must be monotonically increasing.')
427 endif
428 enddo
429 if (value < array(1) .or. value > array(ia)) then
430! if (value < array(1)) frac_index = 1.
431! if (value > array(ia)) frac_index = float(ia)
432 frac_index = -1.0
433 else
434 i=1
435 keep_going = .true.
436 do while (i <= ia .and. keep_going)
437 i = i+1
438 if (value <= array(i)) then
439 frac_index = float(i-1) + (value-array(i-1))/(array(i)-array(i-1))
440 keep_going = .false.
441 endif
442 enddo
443 endif
444 end function frac_index
445
446 !> @brief Return index of nearest point along axis
447 !!
448 !> nearest_index = index of nearest data point within "array" corresponding to
449 !! "value".
450 !!
451 !! inputs:
452 !!
453 !! value = arbitrary data...same units as elements in "array"
454 !! array = array of data points (must be monotonically increasing)
455 !! ia = dimension of "array"
456 !!
457 !! output:
458 !!
459 !! nearest_index = index of nearest data point to "value"
460 !! if "value" is outside the domain of "array" then nearest_index = 1
461 !! or "ia" depending on whether array(1) or array(ia) is
462 !! closest to "value"
463 !!
464 !! note: if "array" is dimensioned array(0:ia) in the calling
465 !! program, then the returned index should be reduced
466 !! by one to account for the zero base.
467 !!
468 !! example:
469 !!
470 !! let model depths be defined by the following:
471 !! parameter (km=5)
472 !! dimension z(km)
473 !! data z /5.0, 10.0, 50.0, 100.0, 250.0/
474 !!
475 !! k1 = nearest_index (12.5, z, km)
476 !! k2 = nearest_index (0.0, z, km)
477 !!
478 !! k1 would be set to 2, and k2 would be set to 1 so that
479 !! z(k1) would be the nearest data point to 12.5 and z(k2) would
480 !! be the nearest data point to 0.0
481 !! @return integer nearest_index
482 function nearest_index (value, array)
483 !=======================================================================
484 !
485 !
486 !=======================================================================
487
488 integer :: nearest_index !< index of nearest data point to "value"
489 !! if "value" is outside the domain of "array" then nearest_index = 1
490 !! or "ia" depending on whether array(1) or array(ia) is
491 !! closest to "value"
492 Integer :: i, ii, unit
493 integer :: ia !< dimension of "array"
494 real :: value !< arbitrary data...same units as elements in "array"
495 real, dimension(:) :: array !< array of data points (must be monotonically increasing)
496 logical keep_going
497
498 ia = size(array(:))
499
500 do i=2,ia
501 if (array(i) < array(i-1)) then
502 unit = stdout()
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:'
507 do ii=1,ia
508 write (unit,*) 'i=',ii, ' array(i)=',array(ii)
509 enddo
510 call mpp_error(fatal,' "nearest_index" array must be monotonically increasing.')
511 endif
512 enddo
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
516 else
517 i=1
518 keep_going = .true.
519 do while (i <= ia .and. keep_going)
520 i = i+1
521 if (value <= array(i)) then
522 nearest_index = i
523 if (array(i)-value > value-array(i-1)) nearest_index = i-1
524 keep_going = .false.
525 endif
526 enddo
527 endif
528 end function nearest_index
529
530 !#############################################################################
531
532 subroutine interp_1d_linear(grid1,grid2,data1,data2)
533
534 real, dimension(:), intent(in) :: grid1, data1, grid2
535 real, dimension(:), intent(inout) :: data2
536
537 integer :: n1, n2, i, n
538 real :: w
539
540 n1 = size(grid1(:))
541 n2 = size(grid2(:))
542
543
544 do i=2,n1
545 if (grid1(i) <= grid1(i-1)) call mpp_error(fatal, 'grid1 not monotonic')
546 enddo
547
548 do i=2,n2
549 if (grid2(i) <= grid2(i-1)) call mpp_error(fatal, 'grid2 not monotonic')
550 enddo
551
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')
554
555 do i=1,n2
556 n = nearest_index(grid2(i),grid1)
557
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)
561 else
562 if(n==1) then
563 data2(i) = data1(n)
564 else
565 w = (grid2(i)-grid1(n-1))/(grid1(n)-grid1(n-1))
566 data2(i) = (1.-w)*data1(n-1) + w*data1(n)
567 endif
568 endif
569 enddo
570
571
572 return
573
574 end subroutine interp_1d_linear
575
576 !###################################################################
577 subroutine interp_1d_cubic_spline(grid1, grid2, data1, data2, yp1, ypn)
578
579 real, dimension(:), intent(in) :: grid1, grid2, data1
580 real, dimension(:), intent(inout) :: data2
581 real, intent(in) :: yp1, ypn
582
583 real, dimension(size(grid1)) :: y2, u
584 real :: sig, p, qn, un, h, a ,b
585 integer :: n, m, i, k, klo, khi
586
587 n = size(grid1(:))
588 m = size(grid2(:))
589
590 do i=2,n
591 if (grid1(i) <= grid1(i-1)) call mpp_error(fatal, 'grid1 not monotonic')
592 enddo
593
594 do i=2,m
595 if (grid2(i) <= grid2(i-1)) call mpp_error(fatal, 'grid2 not monotonic')
596 enddo
597
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')
600
601 if (yp1 >.99e30) then
602 y2(1)=0.
603 u(1)=0.
604 else
605 y2(1)=-0.5
606 u(1)=(3./(grid1(2)-grid1(1)))*((data1(2)-data1(1))/(grid1(2)-grid1(1))-yp1)
607 endif
608
609 do i=2,n-1
610 sig=(grid1(i)-grid1(i-1))/(grid1(i+1)-grid1(i-1))
611 p=sig*y2(i-1)+2.
612 y2(i)=(sig-1.)/p
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
615 enddo
616
617 if (ypn > .99e30) then
618 qn=0.
619 un=0.
620 else
621 qn=0.5
622 un=(3./(grid1(n)-grid1(n-1)))*(ypn-(data1(n)-data1(n-1))/(grid1(n)-grid1(n-1)))
623 endif
624
625 y2(n)=(un-qn*u(n-1))/(qn*y2(n-1)+1.)
626
627 do k=n-1,1,-1
628 y2(k)=y2(k)*y2(k+1)+u(k)
629 enddo
630
631 do k = 1, m
632 n = nearest_index(grid2(k),grid1)
633 if (grid1(n) < grid2(k)) then
634 klo = n
635 else
636 if(n==1) then
637 klo = n
638 else
639 klo = n -1
640 endif
641 endif
642 khi = klo+1
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.
647 enddo
648
649 end subroutine interp_1d_cubic_spline
650
651 !###################################################################
652
653 subroutine interp_1d_1d(grid1,grid2,data1,data2, method, yp1, yp2)
654
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
659
660 real :: y1, y2
661 character(len=32) :: interp_method
662 integer :: k2, ks, ke
663
664 k2 = size(grid2(:))
665
666 interp_method = "linear"
667 if(present(method)) interp_method = method
668 y1 = 1.0e30
669 if(present(yp1)) y1 = yp1
670 y2 = 1.0e30
671 if(present(yp2)) y2 = yp2
672 call find_index(grid1, grid2(1), grid2(k2), ks, ke)
673 select case(trim(interp_method))
674 case("linear")
675 call interp_1d_linear(grid1(ks:ke),grid2,data1(ks:ke),data2)
676 case("cubic_spline")
677 call interp_1d_cubic_spline(grid1(ks:ke),grid2,data1(ks:ke),data2, y1, y2)
678 case default
679 call mpp_error(fatal,"axis_utils: interp_method should be linear or cubic_spline")
680 end select
681
682 return
683
684 end subroutine interp_1d_1d
685
686 !###################################################################
687
688
689 subroutine interp_1d_2d(grid1,grid2,data1,data2)
690
691 real, dimension(:,:), intent(in) :: grid1, data1, grid2
692 real, dimension(:,:), intent(inout) :: data2
693
694 integer :: n1, n2, n, k2, ks, ke
695
696 n1 = size(grid1,1)
697 n2 = size(grid2,1)
698 k2 = size(grid2,2)
699
700 if (n1 /= n2) call mpp_error(fatal,'grid size mismatch')
701
702 do n=1,n1
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,:))
705 enddo
706
707 return
708
709 end subroutine interp_1d_2d
710
711 !###################################################################
712
713 subroutine interp_1d_3d(grid1,grid2,data1,data2, method, yp1, yp2)
714
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
719
720 integer :: n1, n2, m1, m2, k2, n, m
721 real :: y1, y2
722 character(len=32) :: interp_method
723 integer :: ks, ke
724 n1 = size(grid1,1)
725 n2 = size(grid2,1)
726 m1 = size(grid1,2)
727 m2 = size(grid2,2)
728 k2 = size(grid2,3)
729
730 interp_method = "linear"
731 if(present(method)) interp_method = method
732 y1 = 1.0e30
733 if(present(yp1)) y1 = yp1
734 y2 = 1.0e30
735 if(present(yp2)) y2 = yp2
736
737 if (n1 /= n2 .or. m1 /= m2) call mpp_error(fatal,'grid size mismatch')
738
739 select case(trim(interp_method))
740 case("linear")
741 do m=1,m1
742 do n=1,n1
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,:))
745 enddo
746 enddo
747 case("cubic_spline")
748 do m=1,m1
749 do n=1,n1
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)
752 enddo
753 enddo
754 case default
755 call mpp_error(fatal,"axis_utils: interp_method should be linear or cubic_spline")
756 end select
757
758 return
759
760 end subroutine interp_1d_3d
761
762
763 !#####################################################################
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
768
769 integer :: k, nk
770
771 nk = size(grid1(:))
772
773 ks = 0; ke = 0
774 do k = 1, nk-1
775 if(grid1(k) <= xs .and. grid1(k+1) > xs ) then
776 ks = k
777 exit
778 endif
779 enddo
780 do k = nk, 2, -1
781 if(grid1(k) >= xe .and. grid1(k-1) < xe ) then
782 ke = k
783 exit
784 endif
785 enddo
786
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')
789
790 end subroutine find_index
791#endif
792end module axis_utils_mod
793!> @}
794! close documentation grouping
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.
Definition mpp_io.F90:522
Error handler.
Definition mpp.F90:382