47 use mpp_mod,
only: comm_tag_1, comm_tag_2
48 use fms_mod,
only: write_version_number
113 #include<file_version.h> 132 call write_version_number(
"HORIZ_INTERP_CONSERVE_MOD", version)
145 type(horiz_interp_type),
intent(inout) :: Interp
146 real,
intent(in),
dimension(:) :: lon_in , lat_in
147 real,
intent(in),
dimension(:) :: lon_out, lat_out
148 integer,
intent(in),
optional :: verbose
152 real,
dimension(size(lat_out(:))-1,2) :: sph
153 real,
dimension(size(lon_out(:))-1,2) :: theta
154 real,
dimension(size(lat_in(:))) :: slat_in
155 real,
dimension(size(lon_in(:))-1) :: dlon_in
156 real,
dimension(size(lat_in(:))-1) :: dsph_in
157 real,
dimension(size(lon_out(:))-1) :: dlon_out
158 real,
dimension(size(lat_out(:))-1) :: dsph_out
159 real :: blon, fac, hpi, tpi, eps
160 integer :: num_iters = 4
161 integer :: i, j, m, n, nlon_in, nlat_in, nlon_out, nlat_out, &
162 iverbose, m2, n2, iter
164 character(len=64) :: mesg
167 'horiz_interp_conserve_new_1dx1d: horiz_interp_conserve_init is not called')
170 'horiz_interp_conserve_new_1dx1d: great_circle_algorithm is not implemented, contact developer')
172 iverbose = 0;
if (
present(verbose)) iverbose = verbose
180 nlon_in =
size(lon_in(:))-1; nlat_in =
size(lat_in(:))-1
181 nlon_out =
size(lon_out(:))-1; nlat_out =
size(lat_out(:))-1
183 allocate ( interp % facj (nlat_out,2), interp % jlat (nlat_out,2), &
184 interp % faci (nlon_out,2), interp % ilon (nlon_out,2), &
185 interp % area_src (nlon_in, nlat_in), &
186 interp % area_dst (nlon_out, nlat_out) )
192 slat_in(j) = sin(lat_in(j))
196 dsph_in(j) = abs(slat_in(j+1)-slat_in(j))
200 dlon_in(i) = abs(lon_in(i+1)-lon_in(i))
205 if (lat_in(1) > lat_in(nlat_in+1)) s2n = .false.
211 dsph_out(n) = abs(sin(lat_out(n+1))-sin(lat_out(n)))
215 theta(m,1) = lon_out(m)
216 theta(m,2) = lon_out(m+1)
217 dlon_out(m) = abs(lon_out(m+1)-lon_out(m))
220 interp%nlon_src = nlon_in; interp%nlat_src = nlat_in
221 interp%nlon_dst = nlon_out; interp%nlat_dst = nlat_out
228 if (lat_out(n) < lat_out(n+1))
then 229 sph(n,1) = sin(lat_out(n))
230 sph(n,2) = sin(lat_out(n+1))
232 sph(n,1) = sin(lat_out(n+1))
233 sph(n,2) = sin(lat_out(n))
244 if ( (s2n .and. (slat_in(j)-sph(n,n2)) <= eps .and. &
245 (sph(n,n2)-slat_in(j+1)) <= eps) .or. &
246 (.not.s2n .and. (slat_in(j+1)-sph(n,n2)) <= eps .and. &
247 (sph(n,n2)-slat_in(j)) <= eps) )
then 248 interp%jlat(n,n2) = j
250 fac = (sph(n,n2)-slat_in(j))/(slat_in(j+1)-slat_in(j))
252 if (n2 == 1) interp%facj(n,n2) = 1.0 - fac
253 if (n2 == 2) interp%facj(n,n2) = fac
255 if (n2 == 1) interp%facj(n,n2) = fac
256 if (n2 == 2) interp%facj(n,n2) = 1.0 - fac
261 if ( interp%jlat(n,n2) /= 0 )
exit 264 eps = epsilon(sph)*
real(10**iter)
267 if ( interp%jlat(n,n2) == 0 )
then 268 write (mesg,710) n,sph(n,n2)
269 710
format (
': n,sph=',i3,f14.7,40x)
270 call mpp_error(fatal,
'horiz_interp_conserve_mod:no latitude index found'//trim(mesg))
281 if ( blon < lon_in(1) ) blon = blon + tpi
282 if ( blon > lon_in(nlon_in+1) ) blon = blon - tpi
287 if ( (lon_in(i)-blon) <= eps .and. &
288 (blon-lon_in(i+1)) <= eps )
then 289 interp%ilon(m,m2) = i
290 fac = (blon-lon_in(i))/(lon_in(i+1)-lon_in(i))
291 if (m2 == 1) interp%faci(m,m2) = 1.0 - fac
292 if (m2 == 2) interp%faci(m,m2) = fac
296 if ( interp%ilon(m,m2) /= 0 )
exit 299 eps = epsilon(blon)*
real(10**iter)
302 if ( interp%ilon(m,m2) == 0 )
then 303 print *,
'lon_out,blon,blon_in,eps=', &
304 theta(m,m2),blon,lon_in(1),lon_in(nlon_in+1),eps
305 call mpp_error(fatal,
'horiz_interp_conserve_mod: no longitude index found')
314 interp%area_src(i,j) = dlon_in(i) * dsph_in(j)
322 interp%area_dst(m,n) = dlon_out(m) * dsph_out(n)
329 if (iverbose > 2)
then 330 write (*,801) (i,interp%ilon(i,1),interp%ilon(i,2), &
331 interp%faci(i,1),interp%faci(i,2),i=1,nlon_out)
332 write (*,802) (j,interp%jlat(j,1),interp%jlat(j,2), &
333 interp%facj(j,1),interp%facj(j,2),j=1,nlat_out)
334 801
format (/,2x,
'i',4x,
'is',5x,
'ie',4x,
'facis',4x,
'facie', &
336 802
format (/,2x,
'j',4x,
'js',5x,
'je',4x,
'facjs',4x,
'facje', &
346 mask_in, mask_out, verbose)
347 type(horiz_interp_type),
intent(inout) :: Interp
348 real,
intent(in),
dimension(:) :: lon_in , lat_in
349 real,
intent(in),
dimension(:,:) :: lon_out, lat_out
350 real,
intent(in),
optional,
dimension(:,:) :: mask_in
351 real,
intent(inout),
optional,
dimension(:,:) :: mask_out
352 integer,
intent(in),
optional :: verbose
356 integer :: create_xgrid_1DX2D_order1, get_maxxgrid, maxxgrid
357 integer :: create_xgrid_great_circle
358 integer :: nlon_in, nlat_in, nlon_out, nlat_out, nxgrid, i, j
359 real,
dimension(size(lon_in(:))-1, size(lat_in(:))-1) :: mask_src
360 integer,
allocatable,
dimension(:) :: i_src, j_src, i_dst, j_dst
361 real,
allocatable,
dimension(:) :: xgrid_area, clon, clat
362 real,
allocatable,
dimension(:,:) :: dst_area, lon_src, lat_src
363 real,
allocatable,
dimension(:) :: lat_in_flip
364 real,
allocatable,
dimension(:,:) :: mask_src_flip
365 integer :: nincrease, ndecrease
370 'horiz_interp_conserve_new_1dx2d: horiz_interp_conserve_init is not called')
372 if( (
size(lon_out,1) .NE.
size(lat_out,1)) .OR. (
size(lon_out,2) .NE.
size(lat_out,2)) ) &
373 call mpp_error(fatal,
'horiz_interp_conserve_mod: size mismatch between lon_out and lat_out')
374 nlon_in =
size(lon_in(:)) - 1; nlat_in =
size(lat_in(:)) - 1
375 nlon_out =
size(lon_out,1) - 1; nlat_out =
size(lon_out,2) - 1
378 if(
present(mask_in))
then 379 if( (
size(mask_in,1) .NE. nlon_in) .OR. (
size(mask_in,2) .NE. nlat_in))
call mpp_error(fatal, &
380 'horiz_interp_conserve_mod: size mismatch between mask_in and lon_in/lat_in')
384 maxxgrid = get_maxxgrid()
385 allocate(i_src(maxxgrid), j_src(maxxgrid), i_dst(maxxgrid), j_dst(maxxgrid) )
386 allocate( xgrid_area(maxxgrid), dst_area(nlon_out, nlat_out) )
392 if( lat_in(j+1) > lat_in(j) )
then 393 nincrease = nincrease + 1
394 else if ( lat_in(j+1) < lat_in(j) )
then 395 ndecrease = ndecrease + 1
399 if(nincrease == nlat_in)
then 401 else if(ndecrease == nlat_in)
then 404 call mpp_error(fatal,
'horiz_interp_conserve_mod: nlat_in should be equal to nincreaase or ndecrease')
409 allocate(lat_in_flip(nlat_in+1), mask_src_flip(nlon_in,nlat_in))
411 lat_in_flip(j) = lat_in(nlat_in+2-j)
414 mask_src_flip(:,j) = mask_src(:,nlat_in+1-j)
416 nxgrid = create_xgrid_1dx2d_order1(nlon_in, nlat_in, nlon_out, nlat_out, lon_in, lat_in_flip, lon_out, lat_out, &
417 mask_src_flip, i_src, j_src, i_dst, j_dst, xgrid_area)
418 deallocate(lat_in_flip, mask_src_flip)
420 nxgrid = create_xgrid_1dx2d_order1(nlon_in, nlat_in, nlon_out, nlat_out, lon_in, lat_in, lon_out, lat_out, &
421 mask_src, i_src, j_src, i_dst, j_dst, xgrid_area)
424 allocate(lon_src(nlon_in+1,nlat_in+1), lat_src(nlon_in+1,nlat_in+1))
425 allocate(clon(maxxgrid), clat(maxxgrid))
427 allocate(mask_src_flip(nlon_in,nlat_in))
430 lon_src(i,j) = lon_in(i)
431 lat_src(i,j) = lat_in(nlat_in+2-j)
435 mask_src_flip(:,j) = mask_src(:,nlat_in+1-j)
437 nxgrid = create_xgrid_great_circle(nlon_in, nlat_in, nlon_out, nlat_out, lon_src, lat_src, lon_out, lat_out, &
438 mask_src_flip, i_src, j_src, i_dst, j_dst, xgrid_area, clon, clat)
439 deallocate(mask_src_flip)
443 lon_src(i,j) = lon_in(i)
444 lat_src(i,j) = lat_in(j)
447 nxgrid = create_xgrid_great_circle(nlon_in, nlat_in, nlon_out, nlat_out, lon_src, lat_src, lon_out, lat_out, &
448 mask_src, i_src, j_src, i_dst, j_dst, xgrid_area, clon, clat)
450 deallocate(lon_src, lat_src, clon, clat)
452 allocate(interp%i_src(nxgrid), interp%j_src(nxgrid) )
453 allocate(interp%i_dst(nxgrid), interp%j_dst(nxgrid) )
454 allocate(interp%area_frac_dst(nxgrid) )
456 interp%nxgrid = nxgrid
457 interp%i_src = i_src(1:nxgrid)+1
458 interp%j_src = j_src(1:nxgrid)+1
459 if(flip_lat) interp%j_src = nlat_in+1-interp%j_src
460 interp%i_dst = i_dst(1:nxgrid)+1
461 interp%j_dst = j_dst(1:nxgrid)+1
466 dst_area(interp%i_dst(i), interp%j_dst(i)) = dst_area(interp%i_dst(i), interp%j_dst(i)) + xgrid_area(i)
470 interp%area_frac_dst(i) = xgrid_area(i)/dst_area(interp%i_dst(i), interp%j_dst(i) )
472 interp%nlon_src = nlon_in; interp%nlat_src = nlat_in
473 interp%nlon_dst = nlon_out; interp%nlat_dst = nlat_out
474 if(
present(mask_out))
then 475 if( (
size(mask_out,1) .NE. nlon_out) .OR. (
size(mask_out,2) .NE. nlat_out) )
call mpp_error(fatal, &
476 'horiz_interp_conserve_mod: size mismatch between mask_out and lon_out/lat_out')
479 mask_out(interp%i_dst(i),interp%j_dst(i)) = mask_out(interp%i_dst(i),interp%j_dst(i)) + interp%area_frac_dst(i)
483 deallocate(i_src, j_src, i_dst, j_dst, xgrid_area, dst_area )
490 mask_in, mask_out, verbose)
491 type(horiz_interp_type),
intent(inout) :: Interp
492 real,
intent(in),
dimension(:,:) :: lon_in , lat_in
493 real,
intent(in),
dimension(:) :: lon_out, lat_out
494 real,
intent(in),
optional,
dimension(:,:) :: mask_in
495 real,
intent(inout),
optional,
dimension(:,:) :: mask_out
496 integer,
intent(in),
optional :: verbose
500 integer :: create_xgrid_2DX1D_order1, get_maxxgrid, maxxgrid
501 integer :: create_xgrid_great_circle
502 integer :: nlon_in, nlat_in, nlon_out, nlat_out, nxgrid, i, j
503 real,
dimension(size(lon_in,1)-1, size(lon_in,2)-1) :: mask_src
504 integer,
allocatable,
dimension(:) :: i_src, j_src, i_dst, j_dst
505 real,
allocatable,
dimension(:) :: xgrid_area, clon, clat
506 real,
allocatable,
dimension(:,:) :: dst_area, lon_dst, lat_dst
509 'horiz_interp_conserve_new_2dx1d: horiz_interp_conserve_init is not called')
511 if( (
size(lon_in,1) .NE.
size(lat_in,1)) .OR. (
size(lon_in,2) .NE.
size(lat_in,2)) ) &
512 call mpp_error(fatal,
'horiz_interp_conserve_mod: size mismatch between lon_in and lat_in')
513 nlon_in =
size(lon_in,1) - 1; nlat_in =
size(lon_in,2) - 1
514 nlon_out =
size(lon_out(:)) - 1; nlat_out =
size(lat_out(:)) - 1
517 if(
present(mask_in))
then 518 if( (
size(mask_in,1) .NE. nlon_in) .OR. (
size(mask_in,2) .NE. nlat_in))
call mpp_error(fatal, &
519 'horiz_interp_conserve_mod: size mismatch between mask_in and lon_in/lat_in')
523 maxxgrid = get_maxxgrid()
524 allocate(i_src(maxxgrid), j_src(maxxgrid), i_dst(maxxgrid), j_dst(maxxgrid) )
525 allocate( xgrid_area(maxxgrid), dst_area(nlon_out, nlat_out) )
528 nxgrid = create_xgrid_2dx1d_order1(nlon_in, nlat_in, nlon_out, nlat_out, lon_in, lat_in, lon_out, lat_out, &
529 mask_src, i_src, j_src, i_dst, j_dst, xgrid_area)
531 allocate(lon_dst(nlon_out+1, nlat_out+1), lat_dst(nlon_out+1, nlat_out+1) )
532 allocate(clon(maxxgrid), clat(maxxgrid))
535 lon_dst(i,j) = lon_out(i)
536 lat_dst(i,j) = lat_out(j)
539 nxgrid = create_xgrid_great_circle(nlon_in, nlat_in, nlon_out, nlat_out, lon_in, lat_in, lon_dst, lat_dst, &
540 mask_src, i_src, j_src, i_dst, j_dst, xgrid_area, clon, clat)
541 deallocate(lon_dst, lat_dst, clon, clat)
543 allocate(interp%i_src(nxgrid), interp%j_src(nxgrid) )
544 allocate(interp%i_dst(nxgrid), interp%j_dst(nxgrid) )
545 allocate(interp%area_frac_dst(nxgrid) )
547 interp%nxgrid = nxgrid
548 interp%i_src = i_src(1:nxgrid)+1
549 interp%j_src = j_src(1:nxgrid)+1
550 interp%i_dst = i_dst(1:nxgrid)+1
551 interp%j_dst = j_dst(1:nxgrid)+1
556 dst_area(interp%i_dst(i), interp%j_dst(i)) = dst_area(interp%i_dst(i), interp%j_dst(i)) + xgrid_area(i)
560 interp%area_frac_dst(i) = xgrid_area(i)/dst_area(interp%i_dst(i), interp%j_dst(i) )
562 interp%nlon_src = nlon_in; interp%nlat_src = nlat_in
563 interp%nlon_dst = nlon_out; interp%nlat_dst = nlat_out
564 if(
present(mask_out))
then 565 if( (
size(mask_out,1) .NE. nlon_out) .OR. (
size(mask_out,2) .NE. nlat_out) )
call mpp_error(fatal, &
566 'horiz_interp_conserve_mod: size mismatch between mask_out and lon_out/lat_out')
569 mask_out(interp%i_dst(i),interp%j_dst(i)) = mask_out(interp%i_dst(i),interp%j_dst(i)) + interp%area_frac_dst(i)
573 deallocate(i_src, j_src, i_dst, j_dst, xgrid_area, dst_area)
580 mask_in, mask_out, verbose)
581 type(horiz_interp_type),
intent(inout) :: Interp
582 real,
intent(in),
dimension(:,:) :: lon_in , lat_in
583 real,
intent(in),
dimension(:,:) :: lon_out, lat_out
584 real,
intent(in),
optional,
dimension(:,:) :: mask_in
585 real,
intent(inout),
optional,
dimension(:,:) :: mask_out
586 integer,
intent(in),
optional :: verbose
589 integer :: create_xgrid_2DX2D_order1, get_maxxgrid, maxxgrid
590 integer :: create_xgrid_great_circle
591 integer :: nlon_in, nlat_in, nlon_out, nlat_out, nxgrid, i
592 real,
dimension(size(lon_in,1)-1, size(lon_in,2)-1) :: mask_src
593 integer,
allocatable,
dimension(:) :: i_src, j_src, i_dst, j_dst
594 real,
allocatable,
dimension(:) :: xgrid_area, clon, clat
595 real,
allocatable,
dimension(:,:) :: dst_area
598 'horiz_interp_conserve_new_2dx2d: horiz_interp_conserve_init is not called')
600 if( (
size(lon_in,1) .NE.
size(lat_in,1)) .OR. (
size(lon_in,2) .NE.
size(lat_in,2)) ) &
601 call mpp_error(fatal,
'horiz_interp_conserve_mod: size mismatch between lon_in and lat_in')
602 if( (
size(lon_out,1) .NE.
size(lat_out,1)) .OR. (
size(lon_out,2) .NE.
size(lat_out,2)) ) &
603 call mpp_error(fatal,
'horiz_interp_conserve_mod: size mismatch between lon_out and lat_out')
604 nlon_in =
size(lon_in,1) - 1; nlat_in =
size(lon_in,2) - 1
605 nlon_out =
size(lon_out,1) - 1; nlat_out =
size(lon_out,2) - 1
608 if(
present(mask_in))
then 609 if( (
size(mask_in,1) .NE. nlon_in) .OR. (
size(mask_in,2) .NE. nlat_in))
call mpp_error(fatal, &
610 'horiz_interp_conserve_mod: size mismatch between mask_in and lon_in/lat_in')
614 maxxgrid = get_maxxgrid()
615 allocate(i_src(maxxgrid), j_src(maxxgrid), i_dst(maxxgrid), j_dst(maxxgrid) )
616 allocate( xgrid_area(maxxgrid), dst_area(nlon_out, nlat_out) )
619 nxgrid = create_xgrid_2dx2d_order1(nlon_in, nlat_in, nlon_out, nlat_out, lon_in, lat_in, lon_out, lat_out, &
620 mask_src, i_src, j_src, i_dst, j_dst, xgrid_area)
622 allocate(clon(maxxgrid), clat(maxxgrid))
623 nxgrid = create_xgrid_great_circle(nlon_in, nlat_in, nlon_out, nlat_out, lon_in, lat_in, lon_out, lat_out, &
624 mask_src, i_src, j_src, i_dst, j_dst, xgrid_area, clon, clat)
625 deallocate(clon, clat)
628 allocate(interp%i_src(nxgrid), interp%j_src(nxgrid) )
629 allocate(interp%i_dst(nxgrid), interp%j_dst(nxgrid) )
630 allocate(interp%area_frac_dst(nxgrid) )
632 interp%nxgrid = nxgrid
633 interp%i_src = i_src(1:nxgrid)+1
634 interp%j_src = j_src(1:nxgrid)+1
635 interp%i_dst = i_dst(1:nxgrid)+1
636 interp%j_dst = j_dst(1:nxgrid)+1
641 dst_area(interp%i_dst(i), interp%j_dst(i)) = dst_area(interp%i_dst(i), interp%j_dst(i)) + xgrid_area(i)
645 interp%area_frac_dst(i) = xgrid_area(i)/dst_area(interp%i_dst(i), interp%j_dst(i) )
648 interp%nlon_src = nlon_in; interp%nlat_src = nlat_in
649 interp%nlon_dst = nlon_out; interp%nlat_dst = nlat_out
650 if(
present(mask_out))
then 651 if( (
size(mask_out,1) .NE. nlon_out) .OR. (
size(mask_out,2) .NE. nlat_out) )
call mpp_error(fatal, &
652 'horiz_interp_conserve_mod: size mismatch between mask_out and lon_out/lat_out')
655 mask_out(interp%i_dst(i),interp%j_dst(i)) = mask_out(interp%i_dst(i),interp%j_dst(i)) + interp%area_frac_dst(i)
659 deallocate(i_src, j_src, i_dst, j_dst, xgrid_area, dst_area )
710 real,
intent(in),
dimension(:,:) :: data_in
711 real,
intent(out),
dimension(:,:) :: data_out
712 integer,
intent(in),
optional :: verbose
713 real,
intent(in),
dimension(:,:),
optional :: mask_in
714 real,
intent(out),
dimension(:,:),
optional :: mask_out
717 if (
size(data_in,1) /= interp%nlon_src .or.
size(data_in,2) /= interp%nlat_src) &
718 call mpp_error(fatal,
'horiz_interp_conserve_mod: size of input array incorrect')
720 if (
size(data_out,1) /= interp%nlon_dst .or.
size(data_out,2) /= interp%nlat_dst) &
721 call mpp_error(fatal,
'horiz_interp_conserve_mod: size of output array incorrect')
723 select case ( interp%version)
727 if(
present(mask_in) .OR.
present(mask_out) )
call mpp_error(fatal, &
728 'horiz_interp_conserve: for version 2, mask_in and mask_out must be passed in horiz_interp_new, not in horiz_interp')
740 real,
intent(in),
dimension(:,:) :: data_in
741 real,
intent(out),
dimension(:,:) :: data_out
742 integer,
intent(in),
optional :: verbose
743 real,
intent(in),
dimension(:,:),
optional :: mask_in
744 real,
intent(out),
dimension(:,:),
optional :: mask_out
746 integer :: m, n, nlon_in, nlat_in, nlon_out, nlat_out, &
747 miss_in, miss_out, is, ie, js, je, &
749 real :: dsum, wsum, avg_in, min_in, max_in, &
750 avg_out, min_out, max_out, eps, asum, &
751 dwtsum, wtsum, arsum, fis, fie, fjs, fje
753 iverbose = 0;
if (
present(verbose)) iverbose = verbose
757 nlon_in = interp%nlon_src; nlat_in = interp%nlat_src
758 nlon_out = interp%nlon_dst; nlat_out = interp%nlat_dst
760 if (
present(mask_in))
then 761 if ( count(mask_in < -.0001 .or. mask_in > 1.0001) > 0 ) &
762 call mpp_error(fatal,
'horiz_interp_conserve_mod: input mask not between 0,1')
772 if (interp%jlat(n,1) <= interp%jlat(n,2))
then 773 js = interp%jlat(n,1); je = interp%jlat(n,2)
774 fjs = interp%facj(n,1); fje = interp%facj(n,2)
776 js = interp%jlat(n,2); je = interp%jlat(n,1)
777 fjs = interp%facj(n,2); fje = interp%facj(n,1)
782 is = interp%ilon(m,1); ie = interp%ilon(m,2)
783 fis = interp%faci(m,1); fie = interp%faci(m,2)
802 ie = interp%ilon(m,2)
803 fie = interp%faci(m,2)
807 if (
present(mask_in))
then 808 call data_sum ( data_in(is:ie,js:je), interp%area_src(is:ie,js:je), &
809 fis, fie, fjs,fje, dwtsum, wtsum, arsum, mask_in(is:ie,js:je) )
810 else if(
ASSOCIATED(interp%mask_in) )
then 811 call data_sum ( data_in(is:ie,js:je), interp%area_src(is:ie,js:je), &
812 fis, fie, fjs,fje, dwtsum, wtsum, arsum, interp%mask_in(is:ie,js:je) )
814 call data_sum ( data_in(is:ie,js:je), interp%area_src(is:ie,js:je), &
815 fis, fie, fjs,fje, dwtsum, wtsum, arsum )
819 if (wtsum > eps)
then 820 data_out(m,n) = dwtsum/wtsum
821 if (
present(mask_out)) mask_out(m,n) = wtsum/arsum
824 if (
present(mask_out)) mask_out(m,n) = 0.0
834 if (iverbose > 0)
then 838 call stats(data_in, interp%area_src, asum, dsum, wsum, min_in, max_in, miss_in, mask_in)
845 print *,
'horiz_interp stats: input area equals zero ' 848 if (iverbose > 1) print
'(2f16.11)',
'global sum area_in = ', asum, wsum
852 call stats(data_out, interp%area_dst, asum, dsum, wsum, min_out, max_out, miss_out, mask_out)
858 print *,
'horiz_interp stats: output area equals zero ' 861 if (iverbose > 1) print
'(2f16.11)',
'global sum area_out = ', asum, wsum
867 write (*,901) min_in ,max_in ,avg_in
868 if (
present(mask_in))
write (*,903) miss_in
869 write (*,902) min_out,max_out,avg_out
870 if (
present(mask_out))
write (*,903) miss_out
873 900
format (/,1x,10(
'-'),
' output from horiz_interp ',10(
'-'))
874 901
format (
' input: min=',f16.9,
' max=',f16.9,
' avg=',f22.15)
875 902
format (
' output: min=',f16.9,
' max=',f16.9,
' avg=',f22.15)
876 903
format (
' number of missing points = ',i6)
887 real,
intent(in),
dimension(:,:) :: data_in
888 real,
intent(out),
dimension(:,:) :: data_out
889 integer,
intent(in),
optional :: verbose
890 integer :: i, i_src, j_src, i_dst, j_dst
893 do i = 1, interp%nxgrid
894 i_src = interp%i_src(i); j_src = interp%j_src(i)
895 i_dst = interp%i_dst(i); j_dst = interp%j_dst(i)
896 data_out(i_dst, j_dst) = data_out(i_dst, j_dst) + data_in(i_src,j_src)*interp%area_frac_dst(i)
927 select case(interp%version)
929 if(
associated(interp%area_src))
deallocate(interp%area_src)
930 if(
associated(interp%area_dst))
deallocate(interp%area_dst)
931 if(
associated(interp%facj))
deallocate(interp%facj)
932 if(
associated(interp%jlat))
deallocate(interp%jlat)
933 if(
associated(interp%faci))
deallocate(interp%faci)
934 if(
associated(interp%ilon))
deallocate(interp%ilon)
936 if(
associated(interp%i_src))
deallocate(interp%i_src)
937 if(
associated(interp%j_src))
deallocate(interp%j_src)
938 if(
associated(interp%i_dst))
deallocate(interp%i_dst)
939 if(
associated(interp%j_dst))
deallocate(interp%j_dst)
940 if(
associated(interp%area_frac_dst))
deallocate(interp%area_frac_dst)
948 subroutine stats ( dat, area, asum, dsum, wsum, low, high, miss, mask )
949 real,
intent(in) :: dat(:,:), area(:,:)
950 real,
intent(out) :: asum, dsum, wsum, low, high
951 integer,
intent(out) :: miss
952 real,
intent(in),
optional :: mask(:,:)
954 integer :: pe, root_pe, npes, p, buffer_int(1)
955 real :: buffer_real(5)
958 root_pe = mpp_root_pe()
963 if (
present(mask))
then 964 asum = sum(area(:,:))
965 dsum = sum(area(:,:)*dat(:,:)*mask(:,:))
966 wsum = sum(area(:,:)*mask(:,:))
967 miss = count(mask(:,:) <= 0.5)
968 low = minval(dat(:,:),mask=mask(:,:) > 0.5)
969 high = maxval(dat(:,:),mask=mask(:,:) > 0.5)
971 asum = sum(area(:,:))
972 dsum = sum(area(:,:)*dat(:,:))
973 wsum = sum(area(:,:))
975 low = minval(dat(:,:))
976 high = maxval(dat(:,:))
982 if(pe == root_pe)
then 985 call mpp_recv(buffer_real(1),glen=5,from_pe=root_pe+p, tag=comm_tag_1)
986 asum = asum + buffer_real(1)
987 dsum = dsum + buffer_real(2)
988 wsum = wsum + buffer_real(3)
989 low =
min(low, buffer_real(4))
990 high =
max(high, buffer_real(5))
991 call mpp_recv(buffer_int(1),glen=1,from_pe=root_pe+p, tag=comm_tag_2)
992 miss = miss + buffer_int(1)
995 buffer_real(1) = asum
996 buffer_real(2) = dsum
997 buffer_real(3) = wsum
999 buffer_real(5) = high
1001 call mpp_send(buffer_real(1),plen=5,to_pe=root_pe, tag=comm_tag_1)
1002 buffer_int(1) = miss
1003 call mpp_send(buffer_int(1),plen=1,to_pe=root_pe, tag=comm_tag_2)
1006 call mpp_sync_self()
1008 end subroutine stats 1012 subroutine data_sum( data, area, facis, facie, facjs, facje, &
1013 dwtsum, wtsum, arsum, mask )
1017 real,
intent(in),
dimension(:,:) :: data, area
1018 real,
intent(in) :: facis, facie, facjs, facje
1019 real,
intent(inout) :: dwtsum, wtsum, arsum
1020 real,
intent(in),
optional :: mask(:,:)
1028 real,
dimension(size(area,1),size(area,2)) :: wt
1033 id=
size(area,1); jd=
size(area,2)
1036 wt( 1,:)=wt( 1,:)*facis
1037 wt(id,:)=wt(id,:)*facie
1038 wt(:, 1)=wt(:, 1)*facjs
1039 wt(:,jd)=wt(:,jd)*facje
1042 arsum = arsum + asum
1044 if (
present(mask))
then 1046 dwtsum = dwtsum + sum(wt*data)
1047 wtsum = wtsum + sum(wt)
1049 dwtsum = dwtsum + sum(wt*data)
1050 wtsum = wtsum + asum
subroutine horiz_interp_conserve_new_2dx2d(Interp, lon_in, lat_in, lon_out, lat_out, mask_in, mask_out, verbose)
subroutine stats(dat, area, asum, dsum, wsum, low, high, miss, mask)
logical function, public get_great_circle_algorithm()
subroutine horiz_interp_conserve_new_2dx1d(Interp, lon_in, lat_in, lon_out, lat_out, mask_in, mask_out, verbose)
subroutine data_sum(data, area, facis, facie, facjs, facje, dwtsum, wtsum, arsum, mask)
subroutine, public horiz_interp_conserve(Interp, data_in, data_out, verbose, mask_in, mask_out)
subroutine horiz_interp_conserve_new_1dx1d(Interp, lon_in, lat_in, lon_out, lat_out, verbose)
subroutine, public horiz_interp_conserve_del(Interp)
subroutine, public horiz_interp_conserve_init
logical great_circle_algorithm
subroutine horiz_interp_conserve_new_1dx2d(Interp, lon_in, lat_in, lon_out, lat_out, mask_in, mask_out, verbose)
logical module_is_initialized
subroutine horiz_interp_conserve_version2(Interp, data_in, data_out, verbose)
subroutine horiz_interp_conserve_version1(Interp, data_in, data_out, verbose, mask_in, mask_out)
real(fp), parameter, public pi