FV3 Bundle
horiz_interp_bicubic.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  use mpp_mod, only: mpp_error, fatal, stdout, mpp_pe, mpp_root_pe
22  use fms_mod, only: write_version_number
24  use constants_mod, only: pi
25 
26 
27  implicit none
28 
29 ! This module delivers methods for bicubic interpolation from a
30 ! coarse regular grid on a fine regular grid.
31 ! Subroutines
32 !
33 ! bcuint
34 ! bcucof
35 !
36 ! are methods taken from
37 !
38 ! W. H. Press, S. A. Teukolski, W. T. Vetterling and B. P. Flannery,
39 ! Numerical Recipies in FORTRAN, The Art of Scientific Computing.
40 ! Cambridge University Press, 1992
41 !
42 ! written by
43 ! martin.schmidt@io-warnemuende.de (2004)
44 ! revied by
45 ! martin.schmidt@io-warnemuende.de (2004)
46 !
47 ! Version 1.0.0.2005-07-06
48 ! The module is thought to interact with MOM-4.
49 ! Alle benotigten Felder werden extern von MOM verwaltet, da sie
50 ! nicht fur alle interpolierten Daten die gleiche Dimension haben mussen.
51 
52  private
53 
56 
58  module procedure horiz_interp_bicubic_new_1d
59  module procedure horiz_interp_bicubic_new_1d_s
60  end interface
61 
62 ! Include variable "version" to be written to log file.
63 #include<file_version.h>
64  logical :: module_is_initialized = .false.
65  integer :: verbose_bicubic = 0
66 
67 ! Grid variables
68 ! xc, yc : co-ordinates of the coarse grid
69 ! xf, yf : co-ordinates of the fine grid
70 ! fc : variable to be interpolated at the coarse grid
71 ! dfc_x : x-derivative of fc at the coarse grid
72 ! dfc_y : y-derivative of fc at the coarse grid
73 ! dfc_xy : x-y-derivative of fc at the coarse grid
74 ! ff : variable to be interpolated at the fine grid
75 ! dff_x : x-derivative of fc at the fine grid
76 ! dff_y : y-derivative of fc at the fine grid
77 ! dff_xy : x-y-derivative of fc at the fine grid
78 
79 
80  logical :: initialized_bicubic = .false.
81 
82 
83  real, save :: missing = -1e33
84  real :: tpi
85 
86  interface fill_xy
87  module procedure fill_xy
88  end interface
89 
90 
91  contains
92 
93  !#######################################################################
94  ! <SUBROUTINE NAME="horiz_interp_bicubic_init">
95  ! <OVERVIEW>
96  ! writes version number to logfile.out
97  ! </OVERVIEW>
98  ! <DESCRIPTION>
99  ! writes version number to logfile.out
100  ! </DESCRIPTION>
101 
102  subroutine horiz_interp_bicubic_init
104  if(module_is_initialized) return
105  call write_version_number("HORIZ_INTERP_BICUBIC_MOD", version)
106  module_is_initialized = .true.
107  tpi = 2.0*pi
108 
109  end subroutine horiz_interp_bicubic_init
110 
111  ! </SUBROUTINE>
112 
113  !#######################################################################
114  ! <SUBROUTINE NAME="horiz_interp_bicubic_new">
115 
116  ! <OVERVIEW>
117  ! Initialization routine.
118  ! </OVERVIEW>
119  ! <DESCRIPTION>
120  ! Allocates space and initializes a derived-type variable
121  ! that contains pre-computed interpolation indices and weights.
122  ! </DESCRIPTION>
123  ! <TEMPLATE>
124  ! call horiz_interp_bicubic_new ( Interp, lon_in, lat_in, lon_out, lat_out, verbose_bicubic, src_modulo )
125 
126  ! </TEMPLATE>
127  !
128  ! <IN NAME="lon_in" TYPE="real, dimension(:,:)" UNITS="radians">
129  ! Longitude (in radians) for source data grid.
130  ! </IN>
131 
132  ! <IN NAME="lat_in" TYPE="real, dimension(:,:)" UNITS="radians">
133  ! Latitude (in radians) for source data grid.
134  ! </IN>
135 
136  ! <IN NAME="lon_out" TYPE="real, dimension(:,:)" UNITS="radians" >
137  ! Longitude (in radians) for source data grid.
138  ! </IN>
139 
140  ! <IN NAME="lat_out" TYPE="real, dimension(:,:)" UNITS="radians" >
141  ! Latitude (in radians) for source data grid.
142  ! </IN>
143 
144  ! <IN NAME="src_modulo" TYPE="logical, optional">
145  ! logical variable to indicate if the boundary condition along zonal boundary
146  ! is cyclic or not. When true, the zonal boundary condition is cyclic.
147  ! </IN>
148 
149  ! <IN NAME="verbose_bicubic" TYPE="integer, optional" >
150  ! flag for the amount of print output.
151  ! </IN>
152 
153  ! <INOUT NAME="Interp" TYPE="type(horiz_interp_type)" >
154  ! A derived-type variable containing indices and weights used for subsequent
155  ! interpolations. To reinitialize this variable for a different grid-to-grid
156  ! interpolation you must first use the "horiz_interp_bicubic_del" interface.
157  ! </INOUT>
158 
159  subroutine horiz_interp_bicubic_new_1d_s ( Interp, lon_in, lat_in, lon_out, lat_out, &
160  verbose, src_modulo )
162  !-----------------------------------------------------------------------
163  type(horiz_interp_type), intent(inout) :: Interp
164  real, intent(in), dimension(:) :: lon_in , lat_in
165  real, intent(in), dimension(:,:) :: lon_out, lat_out
166  integer, intent(in), optional :: verbose
167  logical, intent(in), optional :: src_modulo
168  integer :: i, j, ip1, im1, jp1, jm1
169  logical :: src_is_modulo
170  integer :: nlon_in, nlat_in, nlon_out, nlat_out
171  integer :: jcl, jcu, icl, icu, jj
172  real :: xz, yz
173  integer :: unit
174 
175  if(present(verbose)) verbose_bicubic = verbose
176  src_is_modulo = .false.
177  if (present(src_modulo)) src_is_modulo = src_modulo
178 
179  if(size(lon_out,1) /= size(lat_out,1) .or. size(lon_out,2) /= size(lat_out,2) ) &
180  call mpp_error(fatal,'horiz_interp_bilinear_mod: when using bilinear ' // &
181  'interplation, the output grids should be geographical grids')
182 
183  !--- get the grid size
184  nlon_in = size(lon_in) ; nlat_in = size(lat_in)
185  nlon_out = size(lon_out,1); nlat_out = size(lat_out,2)
186  interp%nlon_src = nlon_in; interp%nlat_src = nlat_in
187  interp%nlon_dst = nlon_out; interp%nlat_dst = nlat_out
188 ! use wti(:,:,1) for x-derivative, wti(:,:,2) for y-derivative, wti(:,:,3) for xy-derivative
189  allocate ( interp%wti (nlon_in, nlat_in, 3) )
190  allocate ( interp%lon_in (nlon_in) )
191  allocate ( interp%lat_in (nlat_in) )
192  allocate ( interp%rat_x (nlon_out, nlat_out) )
193  allocate ( interp%rat_y (nlon_out, nlat_out) )
194  allocate ( interp%i_lon (nlon_out, nlat_out, 2) )
195  allocate ( interp%j_lat (nlon_out, nlat_out, 2) )
196 
197  interp%lon_in = lon_in
198  interp%lat_in = lat_in
199 
200  if ( verbose_bicubic > 0 ) then
201  unit = stdout()
202  write (unit,'(/,"Initialising bicubic interpolation, interface horiz_interp_bicubic_new_1d_s")')
203  write (unit,'(/," Longitude of coarse grid points (radian): xc(i) i=1, ",i4)') interp%nlon_src
204  write (unit,'(1x,10f10.4)') (interp%lon_in(jj),jj=1,interp%nlon_src)
205  write (unit,'(/," Latitude of coarse grid points (radian): yc(j) j=1, ",i4)') interp%nlat_src
206  write (unit,'(1x,10f10.4)') (interp%lat_in(jj),jj=1,interp%nlat_src)
207  do i=1, interp%nlat_dst
208  write (unit,*)
209  write (unit,'(/," Longitude of fine grid points (radian): xf(i) i=1, ",i4)') interp%nlat_dst
210  write (unit,'(1x,10f10.4)') (lon_out(jj,i),jj=1,interp%nlon_dst)
211  enddo
212  do i=1, interp%nlon_dst
213  write (unit,*)
214  write (unit,'(/," Latitude of fine grid points (radian): yf(j) j=1, ",i4)') interp%nlon_dst
215  write (unit,'(1x,10f10.4)') (lat_out(i,jj),jj=1,interp%nlat_dst)
216  enddo
217  endif
218 
219 
220 !---------------------------------------------------------------------------
221 ! Find the x-derivative. Use central differences and forward or
222 ! backward steps at the boundaries
223 
224  do j=1,nlat_in
225  do i=1,nlon_in
226  ip1=min(i+1,nlon_in)
227  im1=max(i-1,1)
228  interp%wti(i,j,1) = 1./(interp%lon_in(ip1)-interp%lon_in(im1))
229  enddo
230  enddo
231 
232 
233 !---------------------------------------------------------------------------
234 
235 ! Find the y-derivative. Use central differences and forward or
236 ! backward steps at the boundaries
237  do j=1,nlat_in
238  jp1=min(j+1,nlat_in)
239  jm1=max(j-1,1)
240  do i=1,nlon_in
241  interp%wti(i,j,2) = 1./(interp%lat_in(jp1)-interp%lat_in(jm1))
242  enddo
243  enddo
244 
245 !---------------------------------------------------------------------------
246 
247 ! Find the xy-derivative. Use central differences and forward or
248 ! backward steps at the boundaries
249  do j=1,nlat_in
250  jp1=min(j+1,nlat_in)
251  jm1=max(j-1,1)
252  do i=1,nlon_in
253  ip1=min(i+1,nlon_in)
254  im1=max(i-1,1)
255  interp%wti(i,j,3) = 1./((interp%lon_in(ip1)-interp%lon_in(im1))*(interp%lat_in(jp1)-interp%lat_in(jm1)))
256  enddo
257  enddo
258 !---------------------------------------------------------------------------
259 ! Now for each point at the dest-grid find the boundary points of
260 ! the source grid
261  do j=1, nlat_out
262  do i=1,nlon_out
263  yz = lat_out(i,j)
264  xz = lon_out(i,j)
265 
266  jcl = 0
267  jcu = 0
268  if( yz .le. interp%lat_in(1) ) then
269  jcl = 1
270  jcu = 1
271  else if( yz .ge. interp%lat_in(nlat_in) ) then
272  jcl = nlat_in
273  jcu = nlat_in
274  else
275  jcl = indl(interp%lat_in, yz)
276  jcu = indu(interp%lat_in, yz)
277  endif
278 
279  icl = 0
280  icu = 0
281  !--- cyclic condition, do we need to use do while
282  if( xz .gt. interp%lon_in(nlon_in) ) xz = xz - tpi
283  if( xz .le. interp%lon_in(1) ) xz = xz + tpi
284  if( xz .ge. interp%lon_in(nlon_in) ) then
285  icl = nlon_in
286  icu = 1
287  interp%rat_x(i,j) = (xz - interp%lon_in(icl))/(interp%lon_in(icu) - interp%lon_in(icl) + tpi)
288  else
289  icl = indl(interp%lon_in, xz)
290  icu = indu(interp%lon_in, xz)
291  interp%rat_x(i,j) = (xz - interp%lon_in(icl))/(interp%lon_in(icu) - interp%lon_in(icl))
292  endif
293  interp%j_lat(i,j,1) = jcl
294  interp%j_lat(i,j,2) = jcu
295  interp%i_lon(i,j,1) = icl
296  interp%i_lon(i,j,2) = icu
297  if(jcl == jcu) then
298  interp%rat_y(i,j) = 0.0
299  else
300  interp%rat_y(i,j) = (yz - interp%lat_in(jcl))/(interp%lat_in(jcu) - interp%lat_in(jcl))
301  endif
302 ! if(yz.gt.Interp%lat_in(jcu)) call mpp_error(FATAL, ' horiz_interp_bicubic_new_1d_s: yf < ycl, no valid boundary point')
303 ! if(yz.lt.Interp%lat_in(jcl)) call mpp_error(FATAL, ' horiz_interp_bicubic_new_1d_s: yf > ycu, no valid boundary point')
304 ! if(xz.gt.Interp%lon_in(icu)) call mpp_error(FATAL, ' horiz_interp_bicubic_new_1d_s: xf < xcl, no valid boundary point')
305 ! if(xz.lt.Interp%lon_in(icl)) call mpp_error(FATAL, ' horiz_interp_bicubic_new_1d_s: xf > xcu, no valid boundary point')
306  enddo
307  enddo
308  end subroutine horiz_interp_bicubic_new_1d_s
309  ! </SUBROUTINE>
310  subroutine horiz_interp_bicubic_new_1d ( Interp, lon_in, lat_in, lon_out, lat_out, &
311  verbose, src_modulo )
313  !-----------------------------------------------------------------------
314  type(horiz_interp_type), intent(inout) :: Interp
315  real, intent(in), dimension(:) :: lon_in , lat_in
316  real, intent(in), dimension(:) :: lon_out, lat_out
317  integer, intent(in), optional :: verbose
318  logical, intent(in), optional :: src_modulo
319  integer :: i, j, ip1, im1, jp1, jm1
320  logical :: src_is_modulo
321  integer :: nlon_in, nlat_in, nlon_out, nlat_out
322  integer :: jcl, jcu, icl, icu, jj
323  real :: xz, yz
324  integer :: unit
325 
326  if(present(verbose)) verbose_bicubic = verbose
327  src_is_modulo = .false.
328  if (present(src_modulo)) src_is_modulo = src_modulo
329 
330  !--- get the grid size
331  nlon_in = size(lon_in) ; nlat_in = size(lat_in)
332  nlon_out = size(lon_out); nlat_out = size(lat_out)
333  interp%nlon_src = nlon_in; interp%nlat_src = nlat_in
334  interp%nlon_dst = nlon_out; interp%nlat_dst = nlat_out
335  allocate ( interp%wti (nlon_in, nlat_in, 3) )
336  allocate ( interp%lon_in (nlon_in) )
337  allocate ( interp%lat_in (nlat_in) )
338  allocate ( interp%rat_x (nlon_out, nlat_out) )
339  allocate ( interp%rat_y (nlon_out, nlat_out) )
340  allocate ( interp%i_lon (nlon_out, nlat_out, 2) )
341  allocate ( interp%j_lat (nlon_out, nlat_out, 2) )
342 
343  interp%lon_in = lon_in
344  interp%lat_in = lat_in
345 
346  if ( verbose_bicubic > 0 ) then
347  unit = stdout()
348  write (unit,'(/,"Initialising bicubic interpolation, interface horiz_interp_bicubic_new_1d")')
349  write (unit,'(/," Longitude of coarse grid points (radian): xc(i) i=1, ",i4)') interp%nlon_src
350  write (unit,'(1x,10f10.4)') (interp%lon_in(jj),jj=1,interp%nlon_src)
351  write (unit,'(/," Latitude of coarse grid points (radian): yc(j) j=1, ",i4)') interp%nlat_src
352  write (unit,'(1x,10f10.4)') (interp%lat_in(jj),jj=1,interp%nlat_src)
353  write (unit,*)
354  write (unit,'(/," Longitude of fine grid points (radian): xf(i) i=1, ",i4)') interp%nlat_dst
355  write (unit,'(1x,10f10.4)') (lon_out(jj),jj=1,interp%nlon_dst)
356  write (unit,'(/," Latitude of fine grid points (radian): yf(j) j=1, ",i4)') interp%nlon_dst
357  write (unit,'(1x,10f10.4)') (lat_out(jj),jj=1,interp%nlat_dst)
358  endif
359 
360 
361 !---------------------------------------------------------------------------
362 ! Find the x-derivative. Use central differences and forward or
363 ! backward steps at the boundaries
364 
365  do j=1,nlat_in
366  do i=1,nlon_in
367  ip1=min(i+1,nlon_in)
368  im1=max(i-1,1)
369  interp%wti(i,j,1) = 1./(lon_in(ip1)-lon_in(im1))
370  enddo
371  enddo
372 
373 
374 !---------------------------------------------------------------------------
375 
376 ! Find the y-derivative. Use central differences and forward or
377 ! backward steps at the boundaries
378  do j=1,nlat_in
379  jp1=min(j+1,nlat_in)
380  jm1=max(j-1,1)
381  do i=1,nlon_in
382  interp%wti(i,j,2) = 1./(lat_in(jp1)-lat_in(jm1))
383  enddo
384  enddo
385 
386 !---------------------------------------------------------------------------
387 
388 ! Find the xy-derivative. Use central differences and forward or
389 ! backward steps at the boundaries
390  do j=1,nlat_in
391  jp1=min(j+1,nlat_in)
392  jm1=max(j-1,1)
393  do i=1,nlon_in
394  ip1=min(i+1,nlon_in)
395  im1=max(i-1,1)
396  interp%wti(i,j,3) = 1./((lon_in(ip1)-lon_in(im1))*(lat_in(jp1)-lat_in(jm1)))
397  enddo
398  enddo
399 !---------------------------------------------------------------------------
400 ! Now for each point at the dest-grid find the boundary points of
401 ! the source grid
402  do j=1, nlat_out
403  yz = lat_out(j)
404  jcl = 0
405  jcu = 0
406  if( yz .le. lat_in(1) ) then
407  jcl = 1
408  jcu = 1
409  else if( yz .ge. lat_in(nlat_in) ) then
410  jcl = nlat_in
411  jcu = nlat_in
412  else
413  jcl = indl(lat_in, yz)
414  jcu = indu(lat_in, yz)
415  endif
416  do i=1,nlon_out
417  xz = lon_out(i)
418  icl = 0
419  icu = 0
420  !--- cyclic condition, do we need to use do while
421  if( xz .gt. lon_in(nlon_in) ) xz = xz - tpi
422  if( xz .le. lon_in(1) ) xz = xz + tpi
423  if( xz .ge. lon_in(nlon_in) ) then
424  icl = nlon_in
425  icu = 1
426  interp%rat_x(i,j) = (xz - interp%lon_in(icl))/(interp%lon_in(icu) - interp%lon_in(icl) + tpi)
427  else
428  icl = indl(lon_in, xz)
429  icu = indu(lon_in, xz)
430  interp%rat_x(i,j) = (xz - interp%lon_in(icl))/(interp%lon_in(icu) - interp%lon_in(icl))
431  endif
432  icl = indl(lon_in, xz)
433  icu = indu(lon_in, xz)
434  interp%j_lat(i,j,1) = jcl
435  interp%j_lat(i,j,2) = jcu
436  interp%i_lon(i,j,1) = icl
437  interp%i_lon(i,j,2) = icu
438  if(jcl == jcu) then
439  interp%rat_y(i,j) = 0.0
440  else
441  interp%rat_y(i,j) = (yz - interp%lat_in(jcl))/(interp%lat_in(jcu) - interp%lat_in(jcl))
442  endif
443 ! if(yz.gt.lat_in(jcu)) call mpp_error(FATAL, ' horiz_interp_bicubic_new_1d: yf < ycl, no valid boundary point')
444 ! if(yz.lt.lat_in(jcl)) call mpp_error(FATAL, ' horiz_interp_bicubic_new_1d: yf > ycu, no valid boundary point')
445 ! if(xz.gt.lon_in(icu)) call mpp_error(FATAL, ' horiz_interp_bicubic_new_1d: xf < xcl, no valid boundary point')
446 ! if(xz.lt.lon_in(icl)) call mpp_error(FATAL, ' horiz_interp_bicubic_new_1d: xf > xcu, no valid boundary point')
447  enddo
448  enddo
449 
450  end subroutine horiz_interp_bicubic_new_1d
451 
452  subroutine horiz_interp_bicubic( Interp, data_in, data_out, verbose, mask_in, mask_out, missing_value, missing_permit)
453  type(horiz_interp_type), intent(in) :: interp
454  real, intent(in), dimension(:,:) :: data_in
455  real, intent(out), dimension(:,:) :: data_out
456  integer, intent(in), optional :: verbose
457  real, intent(in), dimension(:,:), optional :: mask_in
458  real, intent(out), dimension(:,:), optional :: mask_out
459  real, intent(in), optional :: missing_value
460  integer, intent(in), optional :: missing_permit
461  real :: yz, ycu, ycl
462  real :: xz, xcu, xcl
463  real :: val, val1, val2
464  real, dimension(4) :: y, y1, y2, y12
465  integer :: icl, icu, jcl, jcu
466  integer :: iclp1, icup1, jclp1, jcup1
467  integer :: iclm1, icum1, jclm1, jcum1
468  integer :: i,j
469 
470  if ( present(verbose) ) verbose_bicubic = verbose
471 ! fill_in = .false.
472 ! if ( present(fill) ) fill_in = fill
473 ! use dfc_x and dfc_y as workspace
474 ! if ( fill_in ) call fill_xy(fc(ics:ice,jcs:jce), ics, ice, jcs, jce, maxpass=2)
475 ! where ( data_in .le. missing ) data_in(:,:) = 0.
476 !!
477  do j=1, interp%nlat_dst
478  do i=1, interp%nlon_dst
479  yz = interp%rat_y(i,j)
480  xz = interp%rat_x(i,j)
481  jcl = interp%j_lat(i,j,1)
482  jcu = interp%j_lat(i,j,2)
483  icl = interp%i_lon(i,j,1)
484  icu = interp%i_lon(i,j,2)
485  if( icl > icu ) then
486  iclp1 = icu
487  icum1 = icl
488  xcl = interp%lon_in(icl)
489  xcu = interp%lon_in(icu)+tpi
490  else
491  iclp1 = min(icl+1,interp%nlon_src)
492  icum1 = max(icu-1,1)
493  xcl = interp%lon_in(icl)
494  xcu = interp%lon_in(icu)
495  endif
496  iclm1 = max(icl-1,1)
497  icup1 = min(icu+1,interp%nlon_src)
498  jclp1 = min(jcl+1,interp%nlat_src)
499  jclm1 = max(jcl-1,1)
500  jcup1 = min(jcu+1,interp%nlat_src)
501  jcum1 = max(jcu-1,1)
502  ycl = interp%lat_in(jcl)
503  ycu = interp%lat_in(jcu)
504 ! xcl = Interp%lon_in(icl)
505 ! xcu = Interp%lon_in(icu)
506  y(1) = data_in(icl,jcl)
507  y(2) = data_in(icu,jcl)
508  y(3) = data_in(icu,jcu)
509  y(4) = data_in(icl,jcu)
510  y1(1) = ( data_in(iclp1,jcl) - data_in(iclm1,jcl) ) * interp%wti(icl,jcl,1)
511  y1(2) = ( data_in(icup1,jcl) - data_in(icum1,jcl) ) * interp%wti(icu,jcl,1)
512  y1(3) = ( data_in(icup1,jcu) - data_in(icum1,jcu) ) * interp%wti(icu,jcu,1)
513  y1(4) = ( data_in(iclp1,jcu) - data_in(iclm1,jcu) ) * interp%wti(icl,jcu,1)
514  y2(1) = ( data_in(icl,jclp1) - data_in(icl,jclm1) ) * interp%wti(icl,jcl,2)
515  y2(2) = ( data_in(icu,jclp1) - data_in(icu,jclm1) ) * interp%wti(icu,jcl,2)
516  y2(3) = ( data_in(icu,jcup1) - data_in(icu,jcum1) ) * interp%wti(icu,jcu,2)
517  y2(4) = ( data_in(icl,jcup1) - data_in(icl,jcum1) ) * interp%wti(icl,jcu,2)
518  y12(1)= ( data_in(iclp1,jclp1) + data_in(iclm1,jclm1) - data_in(iclm1,jclp1) &
519  - data_in(iclp1,jclm1) ) * interp%wti(icl,jcl,3)
520  y12(2)= ( data_in(icup1,jclp1) + data_in(icum1,jclm1) - data_in(icum1,jclp1) &
521  - data_in(icup1,jclm1) ) * interp%wti(icu,jcl,3)
522  y12(3)= ( data_in(icup1,jcup1) + data_in(icum1,jcum1) - data_in(icum1,jcup1) &
523  - data_in(icup1,jcum1) ) * interp%wti(icu,jcu,3)
524  y12(4)= ( data_in(iclp1,jcup1) + data_in(iclm1,jcum1) - data_in(iclm1,jcup1) &
525  - data_in(iclp1,jcum1) ) * interp%wti(icl,jcu,3)
526 
527  call bcuint(y,y1,y2,y12,xcl,xcu,ycl,ycu,xz,yz,val,val1,val2)
528  data_out(i,j) = val
529  if(present(mask_out)) mask_out(i,j) = 1.
530 !! dff_x(i,j) = val1
531 !! dff_y(i,j) = val2
532  enddo
533  enddo
534  return
535  end subroutine horiz_interp_bicubic
536 
537 
538 !---------------------------------------------------------------------------
539 
540  subroutine bcuint(y,y1,y2,y12,x1l,x1u,x2l,x2u,t,u,ansy,ansy1,ansy2)
541  real ansy,ansy1,ansy2,x1l,x1u,x2l,x2u,y(4),y1(4),y12(4),y2(4)
542 ! uses bcucof
543  integer i
544  real t,u,c(4,4)
545  call bcucof(y,y1,y2,y12,x1u-x1l,x2u-x2l,c)
546  ansy=0.
547  ansy2=0.
548  ansy1=0.
549  do i=4,1,-1
550  ansy=t*ansy+((c(i,4)*u+c(i,3))*u+c(i,2))*u+c(i,1)
551 ! ansy2=t*ansy2+(3.*c(i,4)*u+2.*c(i,3))*u+c(i,2)
552 ! ansy1=u*ansy1+(3.*c(4,i)*t+2.*c(3,i))*t+c(2,i)
553  enddo
554 ! ansy1=ansy1/(x1u-x1l) ! could be used for accuracy checks
555 ! ansy2=ansy2/(x2u-x2l) ! could be used for accuracy checks
556  return
557 ! (c) copr. 1986-92 numerical recipes software -3#(-)f.
558  end subroutine bcuint
559 !---------------------------------------------------------------------------
560 
561  subroutine bcucof(y,y1,y2,y12,d1,d2,c)
562  real d1,d2,c(4,4),y(4),y1(4),y12(4),y2(4)
563  integer i,j,k,l
564  real d1d2,xx,cl(16),wt(16,16),x(16)
565  save wt
566  data wt/1., 0., -3., 2., 4*0., -3., 0., 9., -6., 2., 0., -6., 4., 8*0., &
567  3., 0., -9., 6., -2., 0., 6., -4., 10*0., 9., -6., 2*0., -6., &
568  4., 2*0., 3., -2., 6*0., -9., 6., 2*0., 6., -4., 4*0., 1., 0., &
569  -3., 2., -2., 0., 6., -4., 1., 0., -3., 2., 8*0., -1., 0., 3., &
570  -2., 1., 0., -3., 2., 10*0., -3., 2., 2*0., 3., -2., 6*0., 3., &
571  -2., 2*0., -6., 4., 2*0., 3., -2., 0., 1., -2., 1., 5*0., -3., &
572  6., -3., 0., 2., -4., 2., 9*0., 3., -6., 3., 0., -2., 4., -2., &
573  10*0., -3., 3., 2*0., 2., -2., 2*0., -1., 1., 6*0., 3., -3., &
574  2*0., -2., 2., 5*0., 1., -2., 1., 0., -2., 4., -2., 0., 1., -2., &
575  1., 9*0., -1., 2., -1., 0., 1., -2., 1., 10*0., 1., -1., 2*0., &
576  -1., 1., 6*0., -1., 1., 2*0., 2., -2., 2*0., -1., 1./
577 
578  d1d2=d1*d2
579  do i=1,4
580  x(i)=y(i)
581  x(i+4)=y1(i)*d1
582  x(i+8)=y2(i)*d2
583  x(i+12)=y12(i)*d1d2
584  enddo
585  do i=1,16
586  xx=0.
587  do k=1,16
588  xx=xx+wt(i,k)*x(k)
589  enddo
590  cl(i)=xx
591  enddo
592  l=0
593  do i=1,4
594  do j=1,4
595  l=l+1
596  c(i,j)=cl(l)
597  enddo
598  enddo
599  return
600 ! (c) copr. 1986-92 numerical recipes software -3#(-)f.
601  end subroutine bcucof
602 
603 !-----------------------------------------------------------------------
604 
605  function indl(xc, xf)
606 ! find the lower neighbour of xf in field xc, return is the index
607  real, intent(in) :: xc(1:)
608  real, intent(in) :: xf
609  integer :: indl
610  integer :: ii
611  indl = 1
612  do ii=1, size(xc)
613  if(xc(ii).gt.xf) return
614  indl = ii
615  enddo
616  call mpp_error(fatal,'Error in indl')
617  return
618  end function indl
619 
620 !-----------------------------------------------------------------------
621 
622  function indu(xc, xf)
623 ! find the upper neighbour of xf in field xc, return is the index
624  real, intent(in) :: xc(1:)
625  real, intent(in) :: xf
626  integer :: indu
627  integer :: ii
628  do ii=1, size(xc)
629  indu = ii
630  if(xc(ii).gt.xf) return
631  enddo
632  call mpp_error(fatal,'Error in indu')
633  return
634  end function indu
635 
636 !-----------------------------------------------------------------------
637 
638  subroutine fill_xy(fi, ics, ice, jcs, jce, mask, maxpass)
639  integer, intent(in) :: ics,ice,jcs,jce
640  real, intent(inout) :: fi(ics:ice,jcs:jce)
641  real, intent(in), optional :: mask(ics:ice,jcs:jce)
642  integer, intent(in) :: maxpass
643  real :: work_old(ics:ice,jcs:jce)
644  real :: work_new(ics:ice,jcs:jce)
645  logical :: ready
646  real :: blank = -1.e30
647  real :: tavr
648  integer :: ipass = 0
649  integer :: inl, inr, jnl, jnu, i, j, is, js, iavr
650 
651 
652  ready = .false.
653 
654  work_new(:,:) = fi(:,:)
655  work_old(:,:) = work_new(:,:)
656  ipass = 0
657  if ( present(mask) ) then
658  do while (.not.ready)
659  ipass = ipass+1
660  ready = .true.
661  do j=jcs, jce
662  do i=ics, ice
663  if (work_old(i,j).le.blank) then
664  tavr=0.
665  iavr=0
666  inl = max(i-1,ics)
667  inr = min(i+1,ice)
668  jnl = max(j-1,jcs)
669  jnu = min(j+1,jce)
670  do js=jnl,jnu
671  do is=inl,inr
672  if (work_old(is,js) .ne. blank .and. mask(is,js).ne.0.) then
673  tavr = tavr + work_old(is,js)
674  iavr = iavr+1
675  endif
676  enddo
677  enddo
678  if (iavr.gt.0) then
679  if (iavr.eq.1) then
680 ! spreading is not allowed if the only valid neighbor is a corner point
681 ! otherwise an ill posed cellular automaton is established leading to
682 ! a spreading of constant values in diagonal direction
683 ! if all corner points are blanked the valid neighbor must be a direct one
684 ! and spreading is allowed
685  if (work_old(inl,jnu).eq.blank.and.&
686  work_old(inr,jnu).eq.blank.and.&
687  work_old(inr,jnl).eq.blank.and.&
688  work_old(inl,jnl).eq.blank) then
689  work_new(i,j)=tavr/iavr
690  ready = .false.
691  endif
692  else
693  work_new(i,j)=tavr/iavr
694  ready = .false.
695  endif
696  endif
697  endif
698  enddo ! j
699  enddo ! i
700 ! save changes made during this pass to work_old
701  work_old(:,:)=work_new(:,:)
702  if(ipass.eq.maxpass) ready=.true.
703  enddo !while (.not.ready)
704  fi(:,:) = work_new(:,:)
705  else
706  do while (.not.ready)
707  ipass = ipass+1
708  ready = .true.
709  do j=jcs, jce
710  do i=ics, ice
711  if (work_old(i,j).le.blank) then
712  tavr=0.
713  iavr=0
714  inl = max(i-1,ics)
715  inr = min(i+1,ice)
716  jnl = max(j-1,jcs)
717  jnu = min(j+1,jce)
718  do is=inl,inr
719  do js=jnl,jnu
720  if (work_old(is,js).gt.blank) then
721  tavr = tavr + work_old(is,js)
722  iavr = iavr+1
723  endif
724  enddo
725  enddo
726  if (iavr.gt.0) then
727  if (iavr.eq.1) then
728 ! spreading is not allowed if the only valid neighbor is a corner point
729 ! otherwise an ill posed cellular automaton is established leading to
730 ! a spreading of constant values in diagonal direction
731 ! if all corner points are blanked the valid neighbor must be a direct one
732 ! and spreading is allowed
733  if (work_old(inl,jnu).le.blank.and. &
734  work_old(inr,jnu).le.blank.and. &
735  work_old(inr,jnl).le.blank.and. &
736  work_old(inl,jnl).le.blank) then
737  work_new(i,j)=tavr/iavr
738  ready = .false.
739  endif
740  else
741  work_new(i,j)=tavr/iavr
742  ready = .false.
743  endif
744  endif
745  endif
746  enddo ! j
747  enddo ! i
748 ! save changes made during this pass to work_old
749  work_old(:,:)=work_new(:,:)
750  if(ipass.eq.maxpass) ready=.true.
751  enddo !while (.not.ready)
752  fi(:,:) = work_new(:,:)
753  endif
754  return
755  end subroutine fill_xy
756 
757  subroutine horiz_interp_bicubic_del( Interp )
759  type(horiz_interp_type), intent(inout) :: interp
760 
761  if(associated(interp%rat_x)) deallocate ( interp%rat_x )
762  if(associated(interp%rat_y)) deallocate ( interp%rat_y )
763  if(associated(interp%lon_in)) deallocate ( interp%lon_in )
764  if(associated(interp%lat_in)) deallocate ( interp%lat_in )
765  if(associated(interp%i_lon)) deallocate ( interp%i_lon )
766  if(associated(interp%j_lat)) deallocate ( interp%j_lat )
767  if(associated(interp%wti)) deallocate ( interp%wti )
768 
769  end subroutine horiz_interp_bicubic_del
770 
771 end module horiz_interp_bicubic_mod
Definition: fms.F90:20
subroutine horiz_interp_bicubic_new_1d(Interp, lon_in, lat_in, lon_out, lat_out, verbose, src_modulo)
subroutine bcucof(y, y1, y2, y12, d1, d2, c)
integer function indu(xc, xf)
Definition: mpp.F90:39
integer function indl(xc, xf)
subroutine, public horiz_interp_bicubic(Interp, data_in, data_out, verbose, mask_in, mask_out, missing_value, missing_permit)
subroutine bcuint(y, y1, y2, y12, x1l, x1u, x2l, x2u, t, u, ansy, ansy1, ansy2)
subroutine, public horiz_interp_bicubic_del(Interp)
subroutine, public horiz_interp_bicubic_init
subroutine horiz_interp_bicubic_new_1d_s(Interp, lon_in, lat_in, lon_out, lat_out, verbose, src_modulo)
#define max(a, b)
Definition: mosaic_util.h:33
#define min(a, b)
Definition: mosaic_util.h:32
real(fp), parameter, public pi