22 use fms_mod,
only: write_version_number
63 #include<file_version.h> 105 call write_version_number(
"HORIZ_INTERP_BICUBIC_MOD", version)
160 verbose, src_modulo )
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
176 src_is_modulo = .false.
177 if (
present(src_modulo)) src_is_modulo = src_modulo
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')
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
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) )
197 interp%lon_in = lon_in
198 interp%lat_in = lat_in
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
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)
212 do i=1, interp%nlon_dst
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)
228 interp%wti(i,j,1) = 1./(interp%lon_in(ip1)-interp%lon_in(im1))
241 interp%wti(i,j,2) = 1./(interp%lat_in(jp1)-interp%lat_in(jm1))
255 interp%wti(i,j,3) = 1./((interp%lon_in(ip1)-interp%lon_in(im1))*(interp%lat_in(jp1)-interp%lat_in(jm1)))
268 if( yz .le. interp%lat_in(1) )
then 271 else if( yz .ge. interp%lat_in(nlat_in) )
then 275 jcl =
indl(interp%lat_in, yz)
276 jcu =
indu(interp%lat_in, yz)
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 287 interp%rat_x(i,j) = (xz - interp%lon_in(icl))/(interp%lon_in(icu) - interp%lon_in(icl) +
tpi)
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))
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
298 interp%rat_y(i,j) = 0.0
300 interp%rat_y(i,j) = (yz - interp%lat_in(jcl))/(interp%lat_in(jcu) - interp%lat_in(jcl))
311 verbose, src_modulo )
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
327 src_is_modulo = .false.
328 if (
present(src_modulo)) src_is_modulo = src_modulo
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) )
343 interp%lon_in = lon_in
344 interp%lat_in = lat_in
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)
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)
369 interp%wti(i,j,1) = 1./(lon_in(ip1)-lon_in(im1))
382 interp%wti(i,j,2) = 1./(lat_in(jp1)-lat_in(jm1))
396 interp%wti(i,j,3) = 1./((lon_in(ip1)-lon_in(im1))*(lat_in(jp1)-lat_in(jm1)))
406 if( yz .le. lat_in(1) )
then 409 else if( yz .ge. lat_in(nlat_in) )
then 413 jcl =
indl(lat_in, yz)
414 jcu =
indu(lat_in, yz)
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 426 interp%rat_x(i,j) = (xz - interp%lon_in(icl))/(interp%lon_in(icu) - interp%lon_in(icl) +
tpi)
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))
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
439 interp%rat_y(i,j) = 0.0
441 interp%rat_y(i,j) = (yz - interp%lat_in(jcl))/(interp%lat_in(jcu) - interp%lat_in(jcl))
452 subroutine horiz_interp_bicubic( Interp, data_in, data_out, verbose, mask_in, mask_out, missing_value, missing_permit)
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
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
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)
488 xcl = interp%lon_in(icl)
489 xcu = interp%lon_in(icu)+
tpi 491 iclp1 =
min(icl+1,interp%nlon_src)
493 xcl = interp%lon_in(icl)
494 xcu = interp%lon_in(icu)
497 icup1 =
min(icu+1,interp%nlon_src)
498 jclp1 =
min(jcl+1,interp%nlat_src)
500 jcup1 =
min(jcu+1,interp%nlat_src)
502 ycl = interp%lat_in(jcl)
503 ycu = interp%lat_in(jcu)
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)
527 call bcuint(y,y1,y2,y12,xcl,xcu,ycl,ycu,xz,yz,val,val1,val2)
529 if(
present(mask_out)) mask_out(i,j) = 1.
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)
545 call bcucof(y,y1,y2,y12,x1u-x1l,x2u-x2l,c)
550 ansy=t*ansy+((c(i,4)*u+c(i,3))*u+c(i,2))*u+c(i,1)
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)
564 real d1d2,xx,cl(16),wt(16,16),x(16)
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./
605 function indl(xc, xf)
607 real,
intent(in) :: xc(1:)
608 real,
intent(in) :: xf
613 if(xc(ii).gt.xf)
return 622 function indu(xc, xf)
624 real,
intent(in) :: xc(1:)
625 real,
intent(in) :: xf
630 if(xc(ii).gt.xf)
return 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)
646 real :: blank = -1.e30
649 integer :: inl, inr, jnl, jnu, i, j, is, js, iavr
654 work_new(:,:) = fi(:,:)
655 work_old(:,:) = work_new(:,:)
657 if (
present(mask) )
then 658 do while (.not.ready)
663 if (work_old(i,j).le.blank)
then 672 if (work_old(is,js) .ne. blank .and. mask(is,js).ne.0.)
then 673 tavr = tavr + work_old(is,js)
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
693 work_new(i,j)=tavr/iavr
701 work_old(:,:)=work_new(:,:)
702 if(ipass.eq.maxpass) ready=.true.
704 fi(:,:) = work_new(:,:)
706 do while (.not.ready)
711 if (work_old(i,j).le.blank)
then 720 if (work_old(is,js).gt.blank)
then 721 tavr = tavr + work_old(is,js)
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
741 work_new(i,j)=tavr/iavr
749 work_old(:,:)=work_new(:,:)
750 if(ipass.eq.maxpass) ready=.true.
752 fi(:,:) = work_new(:,:)
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 )
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)
logical initialized_bicubic
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)
logical module_is_initialized
real(fp), parameter, public pi