FV3 Bundle
horiz_interp_spherical.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 !***********************************************************************
20 
21  ! <CONTACT EMAIL="Matthew.Harrison@noaa.gov"> Matthew Harrison </CONTACT>
22  ! <CONTACT EMAIL="Zhi.Liang@noaa.gov"> Zhi Liang </CONTACT>
23 
24  ! <HISTORY SRC="http://www.gfdl.noaa.gov/fms-cgi-bin/cvsweb.cgi/FMS/"/>
25 
26  ! <OVERVIEW>
27  ! Performs spatial interpolation between grids using inverse-distance-weighted scheme.
28  ! </OVERVIEW>
29 
30  ! <DESCRIPTION>
31  ! This module can interpolate data from rectangular/tripolar grid
32  ! to rectangular/tripolar grid. The interpolation scheme is inverse-distance-weighted
33  ! scheme. There is an optional mask field for missing input data.
34  ! An optional output mask field may be used in conjunction with
35  ! the input mask to show where output data exists.
36  ! </DESCRIPTION>
37 
38  use mpp_mod, only : mpp_error, fatal, warning, stdout
39  use mpp_mod, only : mpp_root_pe, mpp_pe
40  use mpp_mod, only : input_nml_file
41  use fms_mod, only : write_version_number, file_exist, close_file
42  use fms_mod, only : check_nml_error, open_namelist_file
43  use constants_mod, only : pi
45 
46  implicit none
47  private
48 
49 
52 
53  integer, parameter :: max_neighbors = 400
54  real, parameter :: max_dist_default = 0.1 ! radians
55  integer, parameter :: num_nbrs_default = 4
56  real, parameter :: large=1.e20
57  real, parameter :: epsln=1.e-10
58 
59  integer :: pe, root_pe
60 
61 
62  !--- namelist interface
63  !<NAMELIST NAME="horiz_interp_spherical_nml">
64  ! <DATA NAME="search_method" TYPE="character(len=32)">
65  ! indicate the searching method to find the nearest neighbor points. Its value
66  ! can be "radial_search" and "full_search", with default value "radial_search".
67  ! when search_method is "radial_search", the search may be not quite accurate for some cases.
68  ! Normally the search will be ok if you chose suitable max_dist.
69  ! When search_method is "full_search", it will be always accurate, but will be slower
70  ! comparing to "radial_search". Normally these two search algorithm will produce same
71  ! results other than order of operation. "radial_search" are recommended to use.
72  ! The purpose to add "full_search" is in case you think you interpolation results is
73  ! not right, you have other option to verify.
74  ! </DATA>
75  !</NAMELIST>
76 
77  character(len=32) :: search_method = "radial_search" ! or "full_search"
78  namelist /horiz_interp_spherical_nml/ search_method
79 
80  !-----------------------------------------------------------------------
81  ! Include variable "version" to be written to log file.
82 #include<file_version.h>
83  logical :: module_is_initialized = .false.
84 
85 contains
86 
87  !#######################################################################
88  ! <SUBROUTINE NAME="horiz_interp_spherical_init">
89  ! <OVERVIEW>
90  ! writes version number to logfile.out
91  ! </OVERVIEW>
92  ! <DESCRIPTION>
93  ! writes version number to logfile.out
94  ! </DESCRIPTION>
95 
97  integer :: unit, ierr, io
98 
99 
100  if(module_is_initialized) return
101  call write_version_number("horiz_interp_spherical_mod", version)
102 #ifdef INTERNAL_FILE_NML
103  read (input_nml_file, horiz_interp_spherical_nml, iostat=io)
104  ierr = check_nml_error(io,'horiz_interp_spherical_nml')
105 #else
106  if (file_exist('input.nml')) then
107  unit = open_namelist_file( )
108  ierr=1; do while (ierr /= 0)
109  read (unit, nml=horiz_interp_spherical_nml, iostat=io, end=10)
110  ierr = check_nml_error(io,'horiz_interp_spherical_nml') ! also initializes nml error codes
111  enddo
112 10 call close_file (unit)
113  endif
114 #endif
115 
116  module_is_initialized = .true.
117 
118 
119 
120 end subroutine horiz_interp_spherical_init
121 
122  ! </SUBROUTINE>
123 
124  !#######################################################################
125  ! <SUBROUTINE NAME="horiz_interp_spherical_new">
126 
127  ! <OVERVIEW>
128  ! Initialization routine.
129  ! </OVERVIEW>
130  ! <DESCRIPTION>
131  ! Allocates space and initializes a derived-type variable
132  ! that contains pre-computed interpolation indices and weights.
133  ! </DESCRIPTION>
134  ! <TEMPLATE>
135  ! call horiz_interp_spherical_new(Interp, lon_in,lat_in,lon_out,lat_out, num_nbrs, max_dist, src_modulo)
136  ! </TEMPLATE>
137  !
138  ! <IN NAME="lon_in" TYPE="real, dimension(:,:)" UNITS="radians">
139  ! Longitude (in radians) for source data grid.
140  ! </IN>
141 
142  ! <IN NAME="lat_in" TYPE="real, dimension(:,:)" UNITS="radians">
143  ! Latitude (in radians) for source data grid.
144  ! </IN>
145 
146  ! <IN NAME="lon_out" TYPE="real, dimension(:,:)" UNITS="radians" >
147  ! Longitude (in radians) for source data grid.
148  ! </IN>
149 
150  ! <IN NAME="lat_out" TYPE="real, dimension(:,:)" UNITS="radians" >
151  ! Latitude (in radians) for source data grid.
152  ! </IN>
153 
154  ! <IN NAME="num_nbrs" TYPE="integer, optional">
155  ! Number of nearest neighbors for regridding. When number of neighbors within
156  ! the radius max_dist ( namelist variable) is less than num_nbrs, All the neighbors
157  ! will be used to interpolate onto destination grid. when number of neighbors within
158  ! the radius max_dist ( namelist variable) is greater than num_nbrs, at least "num_nbrs"
159  ! neighbors will be used to remap onto destination grid.
160  ! </IN>
161 
162  ! <IN NAME="max_dist" TYPE="real, optional" UNITS="radians">
163  ! Maximum region of influence around destination grid points.
164  ! </IN>
165 
166  ! <IN NAME="src_modulo" TYPE="logical, optional">
167  ! logical variable to indicate if the boundary condition along zonal boundary
168  ! is cyclic or not. When true, the zonal boundary condition is cyclic.
169  ! </IN>
170 
171  ! <INOUT NAME="Interp" TYPE="type(horiz_interp_type)">
172  ! A derived-type variable containing indices and weights used for subsequent
173  ! interpolations. To reinitialize this variable for a different grid-to-grid
174  ! interpolation you must first use the "horiz_interp_del" interface.
175  ! </INOUT>
176 
177  subroutine horiz_interp_spherical_new(Interp, lon_in,lat_in,lon_out,lat_out, &
178  num_nbrs, max_dist, src_modulo)
179  type(horiz_interp_type), intent(inout) :: interp
180  real, intent(in), dimension(:,:) :: lon_in, lat_in, lon_out, lat_out
181  integer, intent(in), optional :: num_nbrs
182  real, optional, intent(in) :: max_dist
183  logical, intent(in), optional :: src_modulo
184 
185  !------local variables ---------------------------------------
186  integer :: i, j, n
187  integer :: map_dst_xsize, map_dst_ysize, map_src_xsize, map_src_ysize
188  integer :: map_src_size, num_neighbors
189  real :: max_src_dist, tpi, hpi
190  logical :: src_is_modulo
191  real :: min_theta_dst, max_theta_dst, min_phi_dst, max_phi_dst
192  real :: min_theta_src, max_theta_src, min_phi_src, max_phi_src
193  integer :: map_src_add(size(lon_out,1),size(lon_out,2),max_neighbors)
194  real :: map_src_dist(size(lon_out,1),size(lon_out,2),max_neighbors)
195  integer :: num_found(size(lon_out,1),size(lon_out,2))
196  integer :: ilon(max_neighbors), jlat(max_neighbors)
197  real, dimension(size(lon_out,1),size(lon_out,2)) :: theta_dst, phi_dst
198  real, dimension(size(lon_in,1)*size(lon_in,2)) :: theta_src, phi_src
199 
200  !--------------------------------------------------------------
201 
202  pe = mpp_pe()
203  root_pe = mpp_root_pe()
204 
205  tpi = 2.0*pi; hpi = 0.5*pi
206 
207  num_neighbors = num_nbrs_default
208  if(present(num_nbrs)) num_neighbors = num_nbrs
209  if (num_neighbors <= 0) call mpp_error(fatal,'horiz_interp_spherical_mod: num_neighbors must be > 0')
210 
211  max_src_dist = max_dist_default
212  if (PRESENT(max_dist)) max_src_dist = max_dist
213  interp%max_src_dist = max_src_dist
214 
215  src_is_modulo = .true.
216  if (PRESENT(src_modulo)) src_is_modulo = src_modulo
217 
218  !--- check the grid size comformable
219  map_dst_xsize=size(lon_out,1);map_dst_ysize=size(lon_out,2)
220  map_src_xsize=size(lon_in,1); map_src_ysize=size(lon_in,2)
221  map_src_size = map_src_xsize*map_src_ysize
222 
223  if (map_dst_xsize /= size(lat_out,1) .or. map_dst_ysize /= size(lat_out,2)) &
224  call mpp_error(fatal,'horiz_interp_spherical_mod: destination grids not conformable')
225  if (map_src_xsize /= size(lat_in,1) .or. map_src_ysize /= size(lat_in,2)) &
226  call mpp_error(fatal,'horiz_interp_spherical_mod: source grids not conformable')
227 
228  theta_src = reshape(lon_in,(/map_src_size/))
229  phi_src = reshape(lat_in,(/map_src_size/))
230  theta_dst(:,:) = lon_out(:,:)
231  phi_dst(:,:) = lat_out(:,:)
232 
233  min_theta_dst=tpi;max_theta_dst=0.0;min_phi_dst=pi;max_phi_dst=-pi
234  min_theta_src=tpi;max_theta_src=0.0;min_phi_src=pi;max_phi_src=-pi
235 
236  where(theta_dst<0.0) theta_dst = theta_dst+tpi
237  where(theta_dst>tpi) theta_dst = theta_dst-tpi
238  where(theta_src<0.0) theta_src = theta_src+tpi
239  where(theta_src>tpi) theta_src = theta_src-tpi
240 
241  where(phi_dst < -hpi) phi_dst = -hpi
242  where(phi_dst > hpi) phi_dst = hpi
243  where(phi_src < -hpi) phi_src = -hpi
244  where(phi_src > hpi) phi_src = hpi
245 
246  do j=1,map_dst_ysize
247  do i=1,map_dst_xsize
248  min_theta_dst = min(min_theta_dst,theta_dst(i,j))
249  max_theta_dst = max(max_theta_dst,theta_dst(i,j))
250  min_phi_dst = min(min_phi_dst,phi_dst(i,j))
251  max_phi_dst = max(max_phi_dst,phi_dst(i,j))
252  enddo
253  enddo
254 
255  do i=1,map_src_size
256  min_theta_src = min(min_theta_src,theta_src(i))
257  max_theta_src = max(max_theta_src,theta_src(i))
258  min_phi_src = min(min_phi_src,phi_src(i))
259  max_phi_src = max(max_phi_src,phi_src(i))
260  enddo
261 
262  if (min_phi_dst < min_phi_src) print *, '=> WARNING: latitute of dest grid exceeds src'
263  if (max_phi_dst > max_phi_src) print *, '=> WARNING: latitute of dest grid exceeds src'
264  ! when src is cyclic, no need to print out the following warning.
265  if(.not. src_is_modulo) then
266  if (min_theta_dst < min_theta_src) print *, '=> WARNING : longitude of dest grid exceeds src'
267  if (max_theta_dst > max_theta_src) print *, '=> WARNING : longitude of dest grid exceeds src'
268  endif
269 
270  ! allocate memory to data type
271  if(ASSOCIATED(interp%i_lon)) then
272  if(size(interp%i_lon,1) .NE. map_dst_xsize .OR. &
273  size(interp%i_lon,2) .NE. map_dst_ysize ) call mpp_error(fatal, &
274  .NE..OR.'horiz_interp_spherical_mod: size(Interp%i_lon(:),1) map_dst_xsize '// &
275  .NE.'size(Interp%i_lon(:),2) map_dst_ysize')
276  else
277  allocate(interp%i_lon(map_dst_xsize,map_dst_ysize,max_neighbors), &
278  interp%j_lat(map_dst_xsize,map_dst_ysize,max_neighbors), &
279  interp%src_dist(map_dst_xsize,map_dst_ysize,max_neighbors), &
280  interp%num_found(map_dst_xsize,map_dst_ysize) )
281  endif
282 
283  map_src_add = 0
284  map_src_dist = large
285  num_found = 0
286 
287  !using radial_search to find the nearest points and corresponding distance.
288 
289  select case(trim(search_method))
290  case ("radial_search") ! will be efficient, but may be not so accurate for some cases
291  call radial_search(theta_src, phi_src, theta_dst, phi_dst, map_src_xsize, map_src_ysize, &
292  map_src_add, map_src_dist, num_found, num_neighbors,max_src_dist,src_is_modulo)
293  case ("full_search") ! always accurate, but less efficient.
294  call full_search(theta_src, phi_src, theta_dst, phi_dst, map_src_add, map_src_dist, &
295  num_found, num_neighbors,max_src_dist )
296  case default
297  call mpp_error(fatal,"horiz_interp_spherical_new: nml search_method = "// &
298  trim(search_method)//" is not a valid namelist option")
299  end select
300 
301  do j=1,map_dst_ysize
302  do i=1,map_dst_xsize
303  do n=1,num_found(i,j)
304  if(map_src_add(i,j,n) == 0) then
305  jlat(n) = 0; ilon(n) = 0
306  else
307  jlat(n) = map_src_add(i,j,n)/map_src_xsize + 1
308  ilon(n) = map_src_add(i,j,n) - (jlat(n)-1)*map_src_xsize
309  if(ilon(n) == 0) then
310  jlat(n) = jlat(n) - 1
311  ilon(n) = map_src_xsize
312  endif
313  endif
314  enddo
315  interp%i_lon(i,j,:) = ilon(:)
316  interp%j_lat(i,j,:) = jlat(:)
317  interp%num_found(i,j) = num_found(i,j)
318  interp%src_dist(i,j,:) = map_src_dist(i,j,:)
319  enddo
320  enddo
321 
322  interp%nlon_src = map_src_xsize; interp%nlat_src = map_src_ysize
323  interp%nlon_dst = map_dst_xsize; interp%nlat_dst = map_dst_ysize
324 
325  return
326 
327  end subroutine horiz_interp_spherical_new
328  ! </SUBROUTINE>
329 
330  !#######################################################################
331  ! <SUBROUTINE NAME="horiz_interp_spherical">
332 
333  ! <OVERVIEW>
334  ! Subroutine for performing the horizontal interpolation between two grids.
335  ! </OVERVIEW>
336  ! <DESCRIPTION>
337  ! Subroutine for performing the horizontal interpolation between two grids.
338  ! horiz_interp_spherical_new must be called before calling this routine.
339  ! </DESCRIPTION>
340  ! <TEMPLATE>
341  ! call horiz_interp_spherical( Interp, data_in, data_out, verbose, mask_in, mask_out, missing_value)
342  ! </TEMPLATE>
343  !
344  ! <IN NAME="Interp" TYPE="type(horiz_interp_type)">
345  ! Derived-type variable containing interpolation indices and weights.
346  ! Returned by a previous call to horiz_interp_spherical_new.
347  ! </IN>
348  ! <IN NAME="data_in" TYPE="real, dimension(:,:)">
349  ! Input data on source grid.
350  ! </IN>
351  ! <IN NAME="verbose" TYPE="integer, optional">
352  ! flag for the amount of print output.
353  ! verbose = 0, no output; = 1, min,max,means; = 2, still more
354  ! </IN>
355  ! <IN NAME="mask_in" TYPE="real, dimension(:,:),optional">
356  ! Input mask, must be the same size as the input data. The real value of
357  ! mask_in must be in the range (0.,1.). Set mask_in=0.0 for data points
358  ! that should not be used or have missing data.
359  ! </IN>
360  ! <IN NAME="missing_value" TYPE="real, optional">
361  ! Use the missing_value to indicate missing data.
362  ! </IN>
363  ! <OUT NAME="data_out" TYPE="real, dimension(:,:)">
364  ! Output data on destination grid.
365  ! </OUT>
366  ! <OUT NAME="mask_out" TYPE="real, dimension(:,:),optional">
367  ! Output mask that specifies whether data was computed.
368  ! </OUT>
369 
370  subroutine horiz_interp_spherical( Interp, data_in, data_out, verbose, mask_in, mask_out, missing_value)
371  type(horiz_interp_type), intent(in) :: interp
372  real, intent(in), dimension(:,:) :: data_in
373  real, intent(out), dimension(:,:) :: data_out
374  integer, intent(in), optional :: verbose
375  real, intent(in), dimension(:,:), optional :: mask_in
376  real, intent(out), dimension(:,:), optional :: mask_out
377  real, intent(in), optional :: missing_value
378 
379  !--- some local variables ----------------------------------------
380  real, dimension(Interp%nlon_dst, Interp%nlat_dst,size(Interp%src_dist,3)) :: wt
381  real, dimension(Interp%nlon_src, Interp%nlat_src) :: mask_src
382  real, dimension(Interp%nlon_dst, Interp%nlat_dst) :: mask_dst
383  integer :: nlon_in, nlat_in, nlon_out, nlat_out, num_found
384  integer :: m, n, i, j, k, miss_in, miss_out, i1, i2, j1, j2, iverbose
385  real :: min_in, max_in, avg_in, min_out, max_out, avg_out, sum
386  !-----------------------------------------------------------------
387 
388  iverbose = 0; if (present(verbose)) iverbose = verbose
389 
390  nlon_in = interp%nlon_src; nlat_in = interp%nlat_src
391  nlon_out = interp%nlon_dst; nlat_out = interp%nlat_dst
392 
393  if(size(data_in,1) .ne. nlon_in .or. size(data_in,2) .ne. nlat_in ) &
394  call mpp_error(fatal,'horiz_interp_spherical_mod: size of input array incorrect')
395 
396  if(size(data_out,1) .ne. nlon_out .or. size(data_out,2) .ne. nlat_out ) &
397  call mpp_error(fatal,'horiz_interp_spherical_mod: size of output array incorrect')
398 
399  mask_src = 1.0; mask_dst = 1.0
400  if(present(mask_in)) mask_src = mask_in
401 
402  do n=1,nlat_out
403  do m=1,nlon_out
404  ! neighbors are sorted nearest to farthest
405  ! check nearest to see if it is a land point
406  num_found = interp%num_found(m,n)
407  if(num_found == 0 ) then
408  mask_dst(m,n) = 0.0
409  else
410  i1 = interp%i_lon(m,n,1); j1 = interp%j_lat(m,n,1)
411  if (mask_src(i1,j1) .lt. 0.5) then
412  mask_dst(m,n) = 0.0
413  endif
414 
415  if(num_found .gt. 1 ) then
416  i2 = interp%i_lon(m,n,2); j2 = interp%j_lat(m,n,2)
417  ! compare first 2 nearest neighbors -- if they are nearly
418  ! equidistant then use this mask for robustness
419  if(abs(interp%src_dist(m,n,2)-interp%src_dist(m,n,1)) .lt. epsln) then
420  if((mask_src(i1,j1) .lt. 0.5)) mask_dst(m,n) = 0.0
421  endif
422  endif
423 
424  sum=0.0
425  do k=1, num_found
426  if(mask_src(interp%i_lon(m,n,k),interp%j_lat(m,n,k)) .lt. 0.5 ) then
427  wt(m,n,k) = 0.0
428  else
429  if (interp%src_dist(m,n,k) <= epsln) then
430  wt(m,n,k) = large
431  sum = sum + large
432  else if(interp%src_dist(m,n,k) <= interp%max_src_dist ) then
433  wt(m,n,k) = 1.0/interp%src_dist(m,n,k)
434  sum = sum+wt(m,n,k)
435  else
436  wt(m,n,k) = 0.0
437  endif
438  endif
439  enddo
440  if (sum > epsln) then
441  do k = 1, num_found
442  wt(m,n,k) = wt(m,n,k)/sum
443  enddo
444  else
445  mask_dst(m,n) = 0.0
446  endif
447  endif
448  enddo
449  enddo
450 
451  data_out = 0.0
452  do n=1,nlat_out
453  do m=1,nlon_out
454  if(mask_dst(m,n) .gt. 0.5) then
455  do k=1, interp%num_found(m,n)
456  i = interp%i_lon(m,n,k)
457  j = interp%j_lat(m,n,k)
458  data_out(m,n) = data_out(m,n)+data_in(i,j)*wt(m,n,k)
459  enddo
460  else
461  if(present(missing_value)) then
462  data_out(m,n) = missing_value
463  else
464  data_out(m,n) = 0.0
465  endif
466  endif
467  enddo
468  enddo
469 
470  if(present(mask_out)) mask_out = mask_dst
471 
472  !***********************************************************************
473  ! compute statistics: minimum, maximum, and mean
474  !-----------------------------------------------------------------------
475 
476  if (iverbose > 0) then
477 
478  ! compute statistics of input data
479 
480  call stats (data_in, min_in, max_in, avg_in, miss_in, missing_value, mask=mask_src)
481 
482  ! compute statistics of output data
483  call stats (data_out, min_out, max_out, avg_out, miss_out, missing_value, mask=mask_dst)
484 
485  !---- output statistics ----
486  ! root_pe have the information of global mean, min and max
487  if(pe == root_pe) then
488  write (*,900)
489  write (*,901) min_in ,max_in, avg_in
490  if (present(mask_in)) write (*,903) miss_in
491  write (*,902) min_out,max_out,avg_out
492  if (present(mask_out)) write (*,903) miss_out
493  endif
494 900 format (/,1x,10('-'),' output from horiz_interp ',10('-'))
495 901 format (' input: min=',f16.9,' max=',f16.9,' avg=',f22.15)
496 902 format (' output: min=',f16.9,' max=',f16.9,' avg=',f22.15)
497 903 format (' number of missing points = ',i6)
498 
499  endif
500 
501  return
502  end subroutine horiz_interp_spherical
503 
504  ! </SUBROUTINE>
505  !#######################################################################
506  subroutine horiz_interp_spherical_wght( Interp, wt, verbose, mask_in, mask_out, missing_value)
507  type(horiz_interp_type), intent(in) :: interp
508  real, intent(out), dimension(:,:,:) :: wt
509  integer, intent(in), optional :: verbose
510  real, intent(in), dimension(:,:), optional :: mask_in
511  real, intent(inout), dimension(:,:), optional :: mask_out
512  real, intent(in), optional :: missing_value
513 
514  !--- some local variables ----------------------------------------
515  real, dimension(Interp%nlon_src, Interp%nlat_src) :: mask_src
516  real, dimension(Interp%nlon_dst, Interp%nlat_dst) :: mask_dst
517  integer :: nlon_in, nlat_in, nlon_out, nlat_out, num_found
518  integer :: m, n, i, j, k, miss_in, miss_out, i1, i2, j1, j2, iverbose
519  real :: min_in, max_in, avg_in, min_out, max_out, avg_out, sum
520  !-----------------------------------------------------------------
521 
522  iverbose = 0; if (present(verbose)) iverbose = verbose
523 
524  nlon_in = interp%nlon_src; nlat_in = interp%nlat_src
525  nlon_out = interp%nlon_dst; nlat_out = interp%nlat_dst
526 
527  mask_src = 1.0; mask_dst = 1.0
528  if(present(mask_in)) mask_src = mask_in
529 
530  do n=1,nlat_out
531  do m=1,nlon_out
532  ! neighbors are sorted nearest to farthest
533  ! check nearest to see if it is a land point
534  num_found = interp%num_found(m,n)
535 
536  if (num_found > num_nbrs_default) then
537  print*,'pe=',mpp_pe(),'num_found=',num_found
538  num_found = num_nbrs_default
539  end if
540 
541  if(num_found == 0 ) then
542  mask_dst(m,n) = 0.0
543  else
544  i1 = interp%i_lon(m,n,1); j1 = interp%j_lat(m,n,1)
545  if (mask_src(i1,j1) .lt. 0.5) then
546  mask_dst(m,n) = 0.0
547  endif
548 
549  if(num_found .gt. 1 ) then
550  i2 = interp%i_lon(m,n,2); j2 = interp%j_lat(m,n,2)
551  ! compare first 2 nearest neighbors -- if they are nearly
552  ! equidistant then use this mask for robustness
553  if(abs(interp%src_dist(m,n,2)-interp%src_dist(m,n,1)) .lt. epsln) then
554  if((mask_src(i1,j1) .lt. 0.5)) mask_dst(m,n) = 0.0
555  endif
556  endif
557 
558  sum=0.0
559  do k=1, num_found
560  if(mask_src(interp%i_lon(m,n,k),interp%j_lat(m,n,k)) .lt. 0.5 ) then
561  wt(m,n,k) = 0.0
562  else
563  if (interp%src_dist(m,n,k) <= epsln) then
564  wt(m,n,k) = large
565  sum = sum + large
566  else if(interp%src_dist(m,n,k) <= interp%max_src_dist ) then
567  wt(m,n,k) = 1.0/interp%src_dist(m,n,k)
568  sum = sum+wt(m,n,k)
569  else
570  wt(m,n,k) = 0.0
571  endif
572  endif
573  enddo
574  if (sum > epsln) then
575  do k = 1, num_found
576  wt(m,n,k) = wt(m,n,k)/sum
577  enddo
578  else
579  mask_dst(m,n) = 0.0
580  endif
581  endif
582  enddo
583  enddo
584 
585  return
586  end subroutine horiz_interp_spherical_wght
587  ! </SUBROUTINE>
588 
589  !#######################################################################
590  ! <SUBROUTINE NAME="horiz_interp_spherical_del">
591 
592  ! <OVERVIEW>
593  ! Deallocates memory used by "horiz_interp_type" variables.
594  ! Must be called before reinitializing with horiz_interp_spherical_new.
595  ! </OVERVIEW>
596  ! <DESCRIPTION>
597  ! Deallocates memory used by "horiz_interp_type" variables.
598  ! Must be called before reinitializing with horiz_interp_spherical_new.
599  ! </DESCRIPTION>
600  ! <TEMPLATE>
601  ! call horiz_interp_spherical_del ( Interp )
602  ! </TEMPLATE>
603 
604  ! <INOUT NAME="Interp" TYPE="horiz_interp_type">
605  ! A derived-type variable returned by previous call
606  ! to horiz_interp_spherical_new. The input variable must have
607  ! allocated arrays. The returned variable will contain
608  ! deallocated arrays.
609  ! </INOUT>
610 
611 
612  subroutine horiz_interp_spherical_del( Interp )
614  type(horiz_interp_type), intent(inout) :: interp
615 
616  if(associated(interp%src_dist)) deallocate(interp%src_dist)
617  if(associated(interp%num_found)) deallocate(interp%num_found)
618  if(associated(interp%i_lon)) deallocate(interp%i_lon)
619  if(associated(interp%j_lat)) deallocate(interp%j_lat)
620 
621  end subroutine horiz_interp_spherical_del
622  ! </SUBROUTINE>
623 
624  !#######################################################################
625 
626 
627  subroutine radial_search(theta_src,phi_src,theta_dst,phi_dst, map_src_xsize, map_src_ysize, &
628  map_src_add, map_src_dist, num_found, num_neighbors,max_src_dist,src_is_modulo)
629  real, intent(in), dimension(:) :: theta_src, phi_src
630  real, intent(in), dimension(:,:) :: theta_dst, phi_dst
631  integer, intent(in) :: map_src_xsize, map_src_ysize
632  integer, intent(out), dimension(:,:,:) :: map_src_add
633  real, intent(out), dimension(:,:,:) :: map_src_dist
634  integer, intent(inout), dimension(:,:) :: num_found
635  integer, intent(in) :: num_neighbors
636  real, intent(in) :: max_src_dist
637  logical, intent(in) :: src_is_modulo
638 
639  !---------- local variables ----------------------------------------
640  integer, parameter :: max_nbrs = 50
641  integer :: i, j, jj, i0, j0, n, l,i_left, i_right
642  integer :: map_dst_xsize, map_dst_ysize
643  integer :: i_left1, i_left2, i_right1, i_right2
644  integer :: map_src_size, step, step_size, bound, bound_start, bound_end
645  logical :: continue_search, result, continue_radial_search
646  real :: d, res
647  !------------------------------------------------------------------
648  map_dst_xsize=size(theta_dst,1);map_dst_ysize=size(theta_dst,2)
649  map_src_size = map_src_xsize*map_src_ysize
650 
651  do j=1,map_dst_ysize
652  do i=1,map_dst_xsize
653  continue_search=.true.
654  step = 1
655  step_size = sqrt(real(map_src_size) )
656  do while (continue_search .and. step_size > 0)
657  do while (step <= map_src_size .and. continue_search)
658  ! count land points as nearest neighbors
659  d = spherical_distance(theta_dst(i,j),phi_dst(i,j),theta_src(step),phi_src(step))
660  if (d <= max_src_dist) then
661  result = update_dest_neighbors(map_src_add(i,j,:),map_src_dist(i,j,:), &
662  step,d, num_found(i,j), num_neighbors )
663  if (result) then
664  n = 0
665  i0 = mod(step,map_src_xsize)
666 
667  if (i0 == 0) i0 = map_src_xsize
668  res = float(step)/float(map_src_xsize)
669  j0 = ceiling(res)
670  continue_radial_search = .true.
671  do while (continue_radial_search)
672  continue_radial_search = .false.
673  n = n+1 ! radial counter
674  if(n > max_nbrs) exit
675  ! ************** left boundary *******************************
676  i_left = i0-n
677  if (i_left <= 0) then
678  if (src_is_modulo) then
679  i_left = map_src_xsize + i_left
680  else
681  i_left = 1
682  endif
683  endif
684 
685  do l = 0, 2*n
686  jj = j0 - n - 1 + l
687  if( jj < 0) then
688  bound = ( 1 - jj )*map_src_xsize - i_left
689  else if ( jj >= map_src_ysize ) then
690  bound = ( 2*map_src_ysize - jj ) * map_src_xsize - i_left
691  else
692  bound = jj * map_src_xsize + i_left
693  endif
694 
695  d = spherical_distance(theta_dst(i,j),phi_dst(i,j),theta_src(bound),phi_src(bound))
696  if(d<=max_src_dist) then
697  result = update_dest_neighbors(map_src_add(i,j,:),map_src_dist(i,j,:), &
698  bound,d, num_found(i,j), num_neighbors)
699  if (result) continue_radial_search = .true.
700  endif
701  enddo
702 
703  ! ***************************right boundary *******************************
704  i_right = i0+n
705  if (i_right > map_src_xsize) then
706  if (src_is_modulo) then
707  i_right = i_right - map_src_xsize
708  else
709  i_right = map_src_xsize
710  endif
711  endif
712 
713  do l = 0, 2*n
714  jj = j0 - n - 1 + l
715  if( jj < 0) then
716  bound = ( 1 - jj )*map_src_xsize - i_right
717  else if ( jj >= map_src_ysize ) then
718  bound = ( 2*map_src_ysize - jj) * map_src_xsize - i_right
719 
720  else
721  bound = jj * map_src_xsize + i_right
722  endif
723 
724  d = spherical_distance(theta_dst(i,j),phi_dst(i,j),theta_src(bound),phi_src(bound))
725  if(d<=max_src_dist) then
726  result = update_dest_neighbors(map_src_add(i,j,:),map_src_dist(i,j,:), &
727  bound,d, num_found(i,j), num_neighbors)
728  if (result) continue_radial_search = .true.
729  endif
730  enddo
731 
732  ! ************************* bottom boundary **********************************
733  i_left2 = 0
734  if( i_left > i_right) then
735  i_left1 = 1
736  i_right1 = i_right
737  i_left2 = i_left
738  i_right2 = map_src_xsize
739  else
740  i_left1 = i_left
741  i_right1 = i_right
742  endif
743 
744  jj = j0 - n - 1
745  if( jj < 0 ) then
746  bound_start = ( 1 - jj)*map_src_xsize - i_right1
747  bound_end = ( 1 - jj)*map_src_xsize - i_left1
748  else
749  bound_start = jj * map_src_xsize + i_left1
750  bound_end = jj * map_src_xsize + i_right1
751  endif
752 
753  bound = bound_start
754  do while (bound <= bound_end)
755  d = spherical_distance(theta_dst(i,j),phi_dst(i,j),theta_src(bound),phi_src(bound))
756  if(d<=max_src_dist) then
757  result = update_dest_neighbors(map_src_add(i,j,:),map_src_dist(i,j,:), &
758  bound,d, num_found(i,j), num_neighbors)
759  if (result) continue_radial_search = .true.
760  endif
761  bound = bound + 1
762 
763  enddo
764 
765  if(i_left2 > 0 ) then
766  if( jj < 0 ) then
767  bound_start = ( 1 - jj)*map_src_xsize - i_right2
768  bound_end = ( 1 - jj)*map_src_xsize - i_left2
769  else
770  bound_start = jj * map_src_xsize + i_left2
771  bound_end = jj * map_src_xsize + i_right2
772  endif
773 
774  bound = bound_start
775  do while (bound <= bound_end)
776  d = spherical_distance(theta_dst(i,j),phi_dst(i,j),theta_src(bound),phi_src(bound))
777  if(d<=max_src_dist) then
778  result = update_dest_neighbors(map_src_add(i,j,:),map_src_dist(i,j,:), &
779  bound,d, num_found(i,j), num_neighbors)
780  if (result) continue_radial_search = .true.
781  endif
782  bound = bound + 1
783  enddo
784  endif
785 
786  ! ************************** top boundary ************************************
787  jj = j0 + n - 1
788  if( jj >= map_src_ysize) then
789  bound_start = ( 2*map_src_ysize - jj ) * map_src_xsize - i_right1
790  bound_end = ( 2*map_src_ysize - jj ) * map_src_xsize - i_left1
791  else
792  bound_start = jj * map_src_xsize + i_left1
793  bound_end = jj * map_src_xsize + i_right1
794  endif
795 
796  bound = bound_start
797  do while (bound <= bound_end)
798  d = spherical_distance(theta_dst(i,j),phi_dst(i,j),theta_src(bound),phi_src(bound))
799  if(d<=max_src_dist) then
800  result = update_dest_neighbors(map_src_add(i,j,:),map_src_dist(i,j,:), &
801  bound,d, num_found(i,j), num_neighbors)
802  if (result) continue_radial_search = .true.
803  endif
804  bound = bound + 1
805  enddo
806 
807  if(i_left2 > 0) then
808  if( jj >= map_src_ysize) then
809  bound_start = ( 2*map_src_ysize - jj ) * map_src_xsize - i_right2
810  bound_end = ( 2*map_src_ysize - jj ) * map_src_xsize - i_left2
811  else
812  bound_start = jj * map_src_xsize + i_left2
813  bound_end = jj * map_src_xsize + i_right2
814  endif
815 
816  bound = bound_start
817  do while (bound <= bound_end)
818  d = spherical_distance(theta_dst(i,j),phi_dst(i,j),theta_src(bound),phi_src(bound))
819  if(d<=max_src_dist) then
820  result = update_dest_neighbors(map_src_add(i,j,:),map_src_dist(i,j,:), &
821  bound,d, num_found(i,j), num_neighbors)
822  if (result) continue_radial_search = .true.
823  endif
824  bound = bound + 1
825  enddo
826  endif
827 
828  enddo
829  continue_search = .false. ! stop looking
830  endif
831  endif
832  step=step+step_size
833  enddo ! search loop
834  step = 1
835  step_size = step_size/2
836  enddo
837  enddo
838  enddo
839 
840  return
841 
842  end subroutine radial_search
843 
844 
845  !#####################################################################
846 
847  function update_dest_neighbors(map_src_add, map_src_dist, src_add,d, num_found, min_nbrs)
849  integer, intent(inout), dimension(:) :: map_src_add
850  real, intent(inout), dimension(:) :: map_src_dist
851  integer, intent(in) :: src_add
852  real, intent(in) :: d
853  integer, intent(inout) :: num_found
854  integer, intent(in) :: min_nbrs
855 
856  logical :: update_dest_neighbors, already_exist = .false.
857 
858  integer :: n,m
859 
860  update_dest_neighbors = .false.
861 
862  n = 0
863  nloop : do while ( n .le. num_found )
864  n = n + 1
865  dist_chk : if (d .le. map_src_dist(n)) then
866  do m=n,num_found
867  if (src_add == map_src_add(m)) then
868  already_exist = .true.
869  exit nloop
870  endif
871  enddo
872  if(num_found < max_neighbors) then
873  num_found = num_found + 1
874  else
875  call mpp_error(fatal,'update_dest_neighbors: '// &
876  'number of neighbor points found is greated than maxium neighbor points' )
877  endif
878  do m=num_found,n+1,-1
879  map_src_add(m) = map_src_add(m-1)
880  map_src_dist(m) = map_src_dist(m-1)
881  enddo
882  map_src_add(n) = src_add
883  map_src_dist(n) = d
884  update_dest_neighbors = .true.
885  if( num_found > min_nbrs ) then
886  if( map_src_dist(num_found) > map_src_dist(num_found-1) ) then
887  num_found = num_found - 1
888  endif
889  if( map_src_dist(min_nbrs+1) > map_src_dist(min_nbrs) ) then
890  num_found = min_nbrs
891  endif
892  endif
893  exit nloop ! n loop
894  endif dist_chk
895  end do nloop
896  if(already_exist) return
897 
898  if( .not. update_dest_neighbors ) then
899  if( num_found < min_nbrs ) then
900  num_found = num_found + 1
901  update_dest_neighbors = .true.
902  map_src_add(num_found) = src_add
903  map_src_dist(num_found) = d
904  endif
905  endif
906 
907 
908  return
909 
910  end function update_dest_neighbors
911 
912  !########################################################################
913 ! function spherical_distance(theta1,phi1,theta2,phi2)
914 
915 ! real, intent(in) :: theta1, phi1, theta2, phi2
916 ! real :: spherical_distance
917 
918 ! real :: r1(3), r2(3), cross(3), s, dot, ang
919 
920  ! this is a simple, enough way to calculate distance on the sphere
921  ! first, construct cartesian vectors r1 and r2
922  ! then calculate the cross-product which is proportional to the area
923  ! between the 2 vectors. The angular distance is arcsin of the
924  ! distancealong the sphere
925  !
926  ! theta is longitude and phi is latitude
927  !
928 
929 
930 ! r1(1) = cos(theta1)*cos(phi1);r1(2)=sin(theta1)*cos(phi1);r1(3)=sin(phi1)
931 ! r2(1) = cos(theta2)*cos(phi2);r2(2)=sin(theta2)*cos(phi2);r2(3)=sin(phi2)
932 
933 ! cross(1) = r1(2)*r2(3)-r1(3)*r2(2)
934 ! cross(2) = r1(3)*r2(1)-r1(1)*r2(3)
935 ! cross(3) = r1(1)*r2(2)-r1(2)*r2(1)
936 
937 ! s = sqrt(cross(1)**2.+cross(2)**2.+cross(3)**2.)
938 
939 ! s = min(s,1.0-epsln)
940 
941 ! dot = r1(1)*r2(1) + r1(2)*r2(2) + r1(3)*r2(3)
942 
943 ! if (dot > 0) then
944 ! ang = asin(s)
945 ! else if (dot < 0) then
946 ! ang = pi + asin(s) !? original is pi - asin(s)
947 ! else
948 ! ang = pi/2.
949 ! endif
950 
951 ! spherical_distance = abs(ang) ! in radians
952 
953 ! return
954 
955 ! end function spherical_distance
956  ! The great cycle distance
957  function spherical_distance(theta1,phi1,theta2,phi2)
959  real, intent(in) :: theta1, phi1, theta2, phi2
960  real :: spherical_distance, dot
961 
962  if(theta1 == theta2 .and. phi1 == phi2) then
963  spherical_distance = 0.0
964  return
965  endif
966 
967  dot = cos(phi1)*cos(phi2)*cos(theta1-theta2) + sin(phi1)*sin(phi2)
968  if(dot > 1. ) dot = 1.
969  if(dot < -1.) dot = -1.
970  spherical_distance = acos(dot)
971 
972  return
973 
974  end function spherical_distance
975 
976 
977  !#######################################################################
978 
979  subroutine full_search(theta_src,phi_src,theta_dst,phi_dst,map_src_add, map_src_dist,num_found, &
980  num_neighbors,max_src_dist)
981  real, intent(in), dimension(:) :: theta_src, phi_src
982  real, intent(in), dimension(:,:) :: theta_dst, phi_dst
983  integer, intent(out), dimension(:,:,:) :: map_src_add
984  real, intent(out), dimension(:,:,:) :: map_src_dist
985  integer, intent(out), dimension(:,:) :: num_found
986  integer, intent(in) :: num_neighbors
987  real, intent(in) :: max_src_dist
988 
989  integer :: i,j,map_src_size, step
990  integer :: map_dst_xsize,map_dst_ysize
991  real :: d
992  logical :: found
993 
994  map_dst_xsize=size(theta_dst,1);map_dst_ysize=size(theta_dst,2)
995  map_src_size =size(theta_src(:))
996 
997  do j=1,map_dst_ysize
998  do i=1,map_dst_xsize
999  do step = 1, map_src_size
1000  d = spherical_distance(theta_dst(i,j),phi_dst(i,j),theta_src(step),phi_src(step))
1001  if( d <= max_src_dist) then
1002  found = update_dest_neighbors(map_src_add(i,j,:),map_src_dist(i,j,:), &
1003  step,d,num_found(i,j), num_neighbors )
1004  endif
1005  enddo
1006  enddo
1007  enddo
1008 
1009  end subroutine full_search
1010 
1011  !#######################################################################
1012 
1013 
1014 end module horiz_interp_spherical_mod
Definition: fms.F90:20
real function spherical_distance(theta1, phi1, theta2, phi2)
integer root_pe
NAMELIST NAME="horiz_interp_spherical_nml">
logical function update_dest_neighbors(map_src_add, map_src_dist, src_add, d, num_found, min_nbrs)
subroutine, public horiz_interp_spherical(Interp, data_in, data_out, verbose, mask_in, mask_out, missing_value)
subroutine, public stats(dat, low, high, avg, miss, missing_value, mask)
subroutine full_search(theta_src, phi_src, theta_dst, phi_dst, map_src_add, map_src_dist, num_found, num_neighbors, max_src_dist)
Definition: mpp.F90:39
integer function, public check_nml_error(IOSTAT, NML_NAME)
Definition: fms.F90:658
subroutine radial_search(theta_src, phi_src, theta_dst, phi_dst, map_src_xsize, map_src_ysize, map_src_add, map_src_dist, num_found, num_neighbors, max_src_dist, src_is_modulo)
character(len=input_str_length), dimension(:), allocatable, target, public input_nml_file
Definition: mpp.F90:1378
subroutine, public horiz_interp_spherical_init
subroutine, public horiz_interp_spherical_new(Interp, lon_in, lat_in, lon_out, lat_out, num_nbrs, max_dist, src_modulo)
subroutine, public horiz_interp_spherical_del(Interp)
************************************************************************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 max(a, b)
Definition: mosaic_util.h:33
subroutine, public horiz_interp_spherical_wght(Interp, wt, verbose, mask_in, mask_out, missing_value)
double dot(const double *p1, const double *p2)
Definition: mosaic_util.c:663
#define min(a, b)
Definition: mosaic_util.h:32
real(fp), parameter, public pi