FV3 Bundle
horiz_interp_conserve.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="Bruce.Wyman@noaa.gov"> Bruce Wyman </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 conservative interpolation
28  ! </OVERVIEW>
29 
30  ! <DESCRIPTION>
31  ! This module can conservatively interpolate data from any logically rectangular grid
32  ! to any rectangular grid. The interpolation scheme is area-averaging
33  ! conservative scheme. There is an optional mask field for missing input data in both
34  ! horiz_interp__conserveinit and horiz_interp_conserve. For efficiency purpose, mask should only be
35  ! kept in horiz_interp_init (will remove the mask in horiz_interp in the future).
36  ! There are 1-D and 2-D version of horiz_interp_conserve_init for 1-D and 2-D grid.
37  ! There is a optional argument mask in horiz_interp_conserve_init_2d and no mask should
38  ! to passed into horiz_interp_conserv. optional argument mask will not be passed into
39  ! horiz_interp_conserve_init_1d and optional argument mask may be passed into
40  ! horiz_interp_conserve (For the purpose of reproduce Memphis??? results).
41  ! An optional output mask field may be used in conjunction with the input mask to show
42  ! where output data exists.
43  ! </DESCRIPTION>
44 
45  use mpp_mod, only: mpp_send, mpp_recv, mpp_pe, mpp_root_pe, mpp_npes
46  use mpp_mod, only: mpp_error, fatal, mpp_sync_self
47  use mpp_mod, only: comm_tag_1, comm_tag_2
48  use fms_mod, only: write_version_number
50  use constants_mod, only: pi
52 
53 
54  implicit none
55  private
56 
57  ! public interface
58 
59 
60  ! <INTERFACE NAME="horiz_interp_conserve_new">
61  ! <OVERVIEW>
62  ! Allocates space and initializes a derived-type variable
63  ! that contains pre-computed interpolation indices and weights.
64  ! </OVERVIEW>
65  ! <DESCRIPTION>
66  ! Allocates space and initializes a derived-type variable
67  ! that contains pre-computed interpolation indices and weights
68  ! for improved performance of multiple interpolations between
69  ! the same grids.
70 
71  ! </DESCRIPTION>
72  ! <IN NAME="lon_in" TYPE="real" DIM="dimension(:), dimension(:,:)" UNITS="radians">
73  ! Longitude (in radians) for source data grid.
74  ! </IN>
75  ! <IN NAME="lat_in" TYPE="real" DIM="dimension(:), dimension(:,:)" UNITS="radians">
76  ! Latitude (in radians) for source data grid.
77  ! </IN>
78  ! <IN NAME="lon_out" TYPE="real" DIM="dimension(:), dimension(:,:)" UNITS="radians" >
79  ! Longitude (in radians) for destination data grid.
80  ! </IN>
81  ! <IN NAME="lat_out" TYPE="real" DIM="dimension(:), dimension(:,:)" UNITS="radians" >
82  ! Latitude (in radians) for destination data grid.
83  ! </IN>
84  ! <IN NAME="verbose" TYPE="integer, optional" >
85  ! flag for the amount of print output.
86  ! </IN>
87  ! <IN NAME="mask_in" TYPE="real, dimension(:,:),optional">
88  ! Input mask. must be the size (size(lon_in)-1, size(lon. The real value of
89  ! mask_in must be in the range (0.,1.). Set mask_in=0.0 for data points
90  ! that should not be used or have missing data.
91  ! </IN>
92  ! <OUT NAME="mask_out" TYPE="real, dimension(:,:),optional">
93  ! Output mask that specifies whether data was computed.
94  ! </OUT>
95  ! <INOUT NAME="Interp" TYPE="type(horiz_interp_type)" >
96  ! A derived-type variable containing indices and weights used for subsequent
97  ! interpolations. To reinitialize this variable for a different grid-to-grid
98  ! interpolation you must first use the "horiz_interp_del" interface.
99  ! </INOUT>
101  module procedure horiz_interp_conserve_new_1dx1d
102  module procedure horiz_interp_conserve_new_1dx2d
103  module procedure horiz_interp_conserve_new_2dx1d
104  module procedure horiz_interp_conserve_new_2dx2d
105  end interface
106  ! </INTERFACE>
109 
110  integer :: pe, root_pe
111  !-----------------------------------------------------------------------
112  ! Include variable "version" to be written to log file.
113 #include<file_version.h>
114  logical :: module_is_initialized = .false.
115 
116  logical :: great_circle_algorithm = .false.
117 
118 contains
119 
120  !#######################################################################
121  ! <SUBROUTINE NAME="horiz_interp_conserve_init">
122  ! <OVERVIEW>
123  ! writes version number to logfile.out
124  ! </OVERVIEW>
125  ! <DESCRIPTION>
126  ! writes version number to logfile.out
127  ! </DESCRIPTION>
128 
129  subroutine horiz_interp_conserve_init
131  if(module_is_initialized) return
132  call write_version_number("HORIZ_INTERP_CONSERVE_MOD", version)
133 
135 
136  module_is_initialized = .true.
137 
138  end subroutine horiz_interp_conserve_init
139 
140  ! </SUBROUTINE>
141 
142  !#######################################################################
143  !<PUBLICROUTINE INTERFACE="horiz_interp_conserve_new">
144  subroutine horiz_interp_conserve_new_1dx1d ( Interp, lon_in, lat_in, lon_out, lat_out, verbose)
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
149 
150  !</PUBLICROUTINE>
151  !-----------------------------------------------------------------------
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
163  logical :: s2n
164  character(len=64) :: mesg
165 
166  if(.not. module_is_initialized) call mpp_error(fatal, &
167  'horiz_interp_conserve_new_1dx1d: horiz_interp_conserve_init is not called')
168 
169  if(great_circle_algorithm) call mpp_error(fatal, &
170  'horiz_interp_conserve_new_1dx1d: great_circle_algorithm is not implemented, contact developer')
171  !-----------------------------------------------------------------------
172  iverbose = 0; if (present(verbose)) iverbose = verbose
173 
174  pe = mpp_pe()
175  root_pe = mpp_root_pe()
176  !-----------------------------------------------------------------------
177  hpi = 0.5*pi
178  tpi = 4.*hpi
179  interp%version = 1
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
182 
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) )
187 
188  !-----------------------------------------------------------------------
189  ! --- set-up for input grid boxes ---
190 
191  do j = 1, nlat_in+1
192  slat_in(j) = sin(lat_in(j))
193  enddo
194 
195  do j = 1, nlat_in
196  dsph_in(j) = abs(slat_in(j+1)-slat_in(j))
197  enddo
198 
199  do i = 1,nlon_in
200  dlon_in(i) = abs(lon_in(i+1)-lon_in(i))
201  enddo
202 
203  ! set south to north flag
204  s2n = .true.
205  if (lat_in(1) > lat_in(nlat_in+1)) s2n = .false.
206 
207  !-----------------------------------------------------------------------
208  ! --- set-up for output grid boxes ---
209 
210  do n = 1, nlat_out
211  dsph_out(n) = abs(sin(lat_out(n+1))-sin(lat_out(n)))
212  enddo
213 
214  do m = 1,nlon_out
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))
218  enddo
219 
220  interp%nlon_src = nlon_in; interp%nlat_src = nlat_in
221  interp%nlon_dst = nlon_out; interp%nlat_dst = nlat_out
222  !***********************************************************************
223 
224  !------ set up latitudinal indexing ------
225  !------ make sure output grid goes south to north ------
226 
227  do n = 1, 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))
231  else
232  sph(n,1) = sin(lat_out(n+1))
233  sph(n,2) = sin(lat_out(n))
234  endif
235  enddo
236 
237  interp%jlat = 0
238  do n2 = 1, 2 ! looping on grid box edges
239  do n = 1, nlat_out ! looping on output latitudes
240  eps = 0.0
241  do iter=1,num_iters
242  ! find indices from input latitudes
243  do j = 1, nlat_in
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
249  ! weight with sin(lat) to exactly conserve area-integral
250  fac = (sph(n,n2)-slat_in(j))/(slat_in(j+1)-slat_in(j))
251  if (s2n) then
252  if (n2 == 1) interp%facj(n,n2) = 1.0 - fac
253  if (n2 == 2) interp%facj(n,n2) = fac
254  else
255  if (n2 == 1) interp%facj(n,n2) = fac
256  if (n2 == 2) interp%facj(n,n2) = 1.0 - fac
257  endif
258  exit
259  endif
260  enddo
261  if ( interp%jlat(n,n2) /= 0 ) exit
262  ! did not find this output grid edge in the input grid
263  ! increase tolerance for multiple passes
264  eps = epsilon(sph)*real(10**iter)
265  enddo
266  ! no match
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))
271  endif
272  enddo
273  enddo
274 
275  !------ set up longitudinal indexing ------
276 
277  interp%ilon = 0
278  do m2 = 1, 2 ! looping on grid box edges
279  do m = 1, nlon_out ! looping on output longitudes
280  blon = theta(m,m2)
281  if ( blon < lon_in(1) ) blon = blon + tpi
282  if ( blon > lon_in(nlon_in+1) ) blon = blon - tpi
283  eps = 0.0
284  do iter=1,num_iters
285  ! find indices from input longitudes
286  do i = 1, nlon_in
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
293  exit
294  endif
295  enddo
296  if ( interp%ilon(m,m2) /= 0 ) exit
297  ! did not find this output grid edge in the input grid
298  ! increase tolerance for multiple passes
299  eps = epsilon(blon)*real(10**iter)
300  enddo
301  ! no match
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')
306  endif
307  enddo
308  enddo
309 
310  ! --- area of input grid boxes ---
311 
312  do j = 1,nlat_in
313  do i = 1,nlon_in
314  interp%area_src(i,j) = dlon_in(i) * dsph_in(j)
315  enddo
316  enddo
317 
318  ! --- area of output grid boxes ---
319 
320  do n = 1, nlat_out
321  do m = 1, nlon_out
322  interp%area_dst(m,n) = dlon_out(m) * dsph_out(n)
323  enddo
324  enddo
325 
326  !-----------------------------------------------------------------------
327  ! this output may be quite lengthy and is not recommended
328  ! when using more than one processor
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', &
335  /,(i4,2i7,2f10.5))
336 802 format (/,2x,'j',4x,'js',5x,'je',4x,'facjs',4x,'facje', &
337  /,(i4,2i7,2f10.5))
338  endif
339  !-----------------------------------------------------------------------
340 
341  end subroutine horiz_interp_conserve_new_1dx1d
342 
343  !#######################################################################
344  !<PUBLICROUTINE INTERFACE="horiz_interp_conserve_new">
345  subroutine horiz_interp_conserve_new_1dx2d ( Interp, lon_in, lat_in, lon_out, lat_out, &
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
353 
354  !</PUBLICROUTINE>
355 
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
366  logical :: flip_lat
367 
368 
369  if(.not. module_is_initialized) call mpp_error(fatal, &
370  'horiz_interp_conserve_new_1dx2d: horiz_interp_conserve_init is not called')
371 
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
376 
377  mask_src = 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')
381  mask_src = mask_in
382  end if
383 
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) )
387 
388  !--- check if source latitude is flipped
389  nincrease = 0
390  ndecrease = 0
391  do j = 1, nlat_in
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
396  endif
397  enddo
398 
399  if(nincrease == nlat_in) then
400  flip_lat = .false.
401  else if(ndecrease == nlat_in) then
402  flip_lat = .true.
403  else
404  call mpp_error(fatal, 'horiz_interp_conserve_mod: nlat_in should be equal to nincreaase or ndecrease')
405  endif
406 
407  if( .not. great_circle_algorithm ) then
408  if(flip_lat) then
409  allocate(lat_in_flip(nlat_in+1), mask_src_flip(nlon_in,nlat_in))
410  do j = 1, nlat_in+1
411  lat_in_flip(j) = lat_in(nlat_in+2-j)
412  enddo
413  do j = 1, nlat_in
414  mask_src_flip(:,j) = mask_src(:,nlat_in+1-j)
415  enddo
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)
419  else
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)
422  endif
423  else
424  allocate(lon_src(nlon_in+1,nlat_in+1), lat_src(nlon_in+1,nlat_in+1))
425  allocate(clon(maxxgrid), clat(maxxgrid))
426  if(flip_lat) then
427  allocate(mask_src_flip(nlon_in,nlat_in))
428  do j = 1, nlat_in+1
429  do i = 1, nlon_in+1
430  lon_src(i,j) = lon_in(i)
431  lat_src(i,j) = lat_in(nlat_in+2-j)
432  enddo
433  enddo
434  do j = 1, nlat_in
435  mask_src_flip(:,j) = mask_src(:,nlat_in+1-j)
436  enddo
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)
440  else
441  do j = 1, nlat_in+1
442  do i = 1, nlon_in+1
443  lon_src(i,j) = lon_in(i)
444  lat_src(i,j) = lat_in(j)
445  enddo
446  enddo
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)
449  endif
450  deallocate(lon_src, lat_src, clon, clat)
451  endif
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) )
455  interp%version = 2
456  interp%nxgrid = nxgrid
457  interp%i_src = i_src(1:nxgrid)+1 ! in C, the starting index is 0
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
462 
463  ! sum over exchange grid area to get destination grid area
464  dst_area = 0.
465  do i = 1, nxgrid
466  dst_area(interp%i_dst(i), interp%j_dst(i)) = dst_area(interp%i_dst(i), interp%j_dst(i)) + xgrid_area(i)
467  end do
468 
469  do i = 1, nxgrid
470  interp%area_frac_dst(i) = xgrid_area(i)/dst_area(interp%i_dst(i), interp%j_dst(i) )
471  end do
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')
477  mask_out = 0.0
478  do i = 1, nxgrid
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)
480  end do
481  end if
482 
483  deallocate(i_src, j_src, i_dst, j_dst, xgrid_area, dst_area )
484 
485  end subroutine horiz_interp_conserve_new_1dx2d
486 
487  !#######################################################################
488  !<PUBLICROUTINE INTERFACE="horiz_interp_conserve_new">
489  subroutine horiz_interp_conserve_new_2dx1d ( Interp, lon_in, lat_in, lon_out, lat_out, &
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
497 
498  !</PUBLICROUTINE>
499 
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
507 
508  if(.not. module_is_initialized) call mpp_error(fatal, &
509  'horiz_interp_conserve_new_2dx1d: horiz_interp_conserve_init is not called')
510 
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
515 
516  mask_src = 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')
520  mask_src = mask_in
521  end if
522 
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) )
526 
527  if( .not. great_circle_algorithm ) then
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)
530  else
531  allocate(lon_dst(nlon_out+1, nlat_out+1), lat_dst(nlon_out+1, nlat_out+1) )
532  allocate(clon(maxxgrid), clat(maxxgrid))
533  do j = 1, nlat_out+1
534  do i = 1, nlon_out+1
535  lon_dst(i,j) = lon_out(i)
536  lat_dst(i,j) = lat_out(j)
537  enddo
538  enddo
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)
542  endif
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) )
546  interp%version = 2
547  interp%nxgrid = nxgrid
548  interp%i_src = i_src(1:nxgrid)+1 ! in C, the starting index is 0
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
552 
553  ! sum over exchange grid area to get destination grid area
554  dst_area = 0.
555  do i = 1, nxgrid
556  dst_area(interp%i_dst(i), interp%j_dst(i)) = dst_area(interp%i_dst(i), interp%j_dst(i)) + xgrid_area(i)
557  end do
558 
559  do i = 1, nxgrid
560  interp%area_frac_dst(i) = xgrid_area(i)/dst_area(interp%i_dst(i), interp%j_dst(i) )
561  end do
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')
567  mask_out = 0.0
568  do i = 1, nxgrid
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)
570  end do
571  end if
572 
573  deallocate(i_src, j_src, i_dst, j_dst, xgrid_area, dst_area)
574 
575  end subroutine horiz_interp_conserve_new_2dx1d
576 
577  !#######################################################################
578  !<PUBLICROUTINE INTERFACE="horiz_interp_conserve_new">
579  subroutine horiz_interp_conserve_new_2dx2d ( Interp, lon_in, lat_in, lon_out, lat_out, &
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
587  !</PUBLICROUTINE>
588 
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
596 
597  if(.not. module_is_initialized) call mpp_error(fatal, &
598  'horiz_interp_conserve_new_2dx2d: horiz_interp_conserve_init is not called')
599 
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
606 
607  mask_src = 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')
611  mask_src = mask_in
612  end if
613 
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) )
617 
618  if( .not. great_circle_algorithm ) then
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)
621  else
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)
626  endif
627 
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) )
631  interp%version = 2
632  interp%nxgrid = nxgrid
633  interp%i_src = i_src(1:nxgrid)+1 ! in C, the starting index is 0
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
637 
638  ! sum over exchange grid area to get destination grid area
639  dst_area = 0.
640  do i = 1, nxgrid
641  dst_area(interp%i_dst(i), interp%j_dst(i)) = dst_area(interp%i_dst(i), interp%j_dst(i)) + xgrid_area(i)
642  end do
643 
644  do i = 1, nxgrid
645  interp%area_frac_dst(i) = xgrid_area(i)/dst_area(interp%i_dst(i), interp%j_dst(i) )
646  end do
647 
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')
653  mask_out = 0.0
654  do i = 1, nxgrid
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)
656  end do
657  end if
658 
659  deallocate(i_src, j_src, i_dst, j_dst, xgrid_area, dst_area )
660 
661  end subroutine horiz_interp_conserve_new_2dx2d
662 
663  !########################################################################
664  ! <SUBROUTINE NAME="horiz_interp_conserve">
665 
666  ! <OVERVIEW>
667  ! Subroutine for performing the horizontal interpolation between two grids.
668  ! </OVERVIEW>
669  ! <DESCRIPTION>
670  ! Subroutine for performing the horizontal interpolation between two grids.
671  ! horiz_interp_conserve_new must be called before calling this routine.
672  ! </DESCRIPTION>
673  ! <TEMPLATE>
674  ! call horiz_interp_conserve ( Interp, data_in, data_out, verbose, mask_in, mask_out)
675  ! </TEMPLATE>
676  !
677  ! <IN NAME="Interp" TYPE="type(horiz_interp_type)">
678  ! Derived-type variable containing interpolation indices and weights.
679  ! Returned by a previous call to horiz_interp_new.
680  ! </IN>
681  ! <IN NAME="data_in" TYPE="real, dimension(:,:)">
682  ! Input data on source grid.
683  ! </IN>
684 
685  ! <IN NAME="verbose" TYPE="integer, optional">
686  ! flag for the amount of print output.
687  ! verbose = 0, no output; = 1, min,max,means; = 2, still more
688  ! </IN>
689  ! <IN NAME="mask_in" TYPE="real, dimension(:,:),optional">
690  ! Input mask, must be the same size as the input data. The real value of
691  ! mask_in must be in the range (0.,1.). Set mask_in=0.0 for data points
692  ! that should not be used or have missing data. mask_in will be applied only
693  ! when horiz_interp_conserve_new_1d is called. mask_in will be passed into
694  ! horiz_interp_conserve_new_2d.
695  ! </IN>
696 
697  ! <OUT NAME="data_out" TYPE="real, dimension(:,:)">
698  ! Output data on destination grid.
699  ! </OUT>
700  ! <OUT NAME="mask_out" TYPE="real, dimension(:,:),optional">
701  ! Output mask that specifies whether data was computed. mask_out will be computed only
702  ! when horiz_interp_conserve_new_1d is called. mask_out will be computed in
703  ! horiz_interp_conserve_new_2d.
704  ! </OUT>
705 
706  subroutine horiz_interp_conserve ( Interp, data_in, data_out, verbose, &
707  mask_in, mask_out)
708  !-----------------------------------------------------------------------
709  type(horiz_interp_type), intent(in) :: interp
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
715 
716  ! --- error checking ---
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')
719 
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')
722 
723  select case ( interp%version)
724  case (1)
725  call horiz_interp_conserve_version1(interp, data_in, data_out, verbose, mask_in, mask_out)
726  case (2)
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')
729  call horiz_interp_conserve_version2(interp, data_in, data_out, verbose)
730  end select
731 
732  end subroutine horiz_interp_conserve
733  ! </SUBROUTINE>
734 
735  !##############################################################################
736  subroutine horiz_interp_conserve_version1 ( Interp, data_in, data_out, verbose, &
737  mask_in, mask_out)
738  !-----------------------------------------------------------------------
739  type(horiz_interp_type), intent(in) :: 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
745  !----------local variables----------------------------------------------------
746  integer :: m, n, nlon_in, nlat_in, nlon_out, nlat_out, &
747  miss_in, miss_out, is, ie, js, je, &
748  np, npass, iverbose
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
752  !-----------------------------------------------------------------------
753  iverbose = 0; if (present(verbose)) iverbose = verbose
754 
755  eps = epsilon(wtsum)
756 
757  nlon_in = interp%nlon_src; nlat_in = interp%nlat_src
758  nlon_out = interp%nlon_dst; nlat_out = interp%nlat_dst
759 
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')
763  endif
764 
765  !-----------------------------------------------------------------------
766  !---- loop through output grid boxes ----
767 
768  data_out = 0.0
769  do n = 1, nlat_out
770  ! latitude window
771  ! setup ascending latitude indices and weights
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)
775  else
776  js = interp%jlat(n,2); je = interp%jlat(n,1)
777  fjs = interp%facj(n,2); fje = interp%facj(n,1)
778  endif
779 
780  do m = 1, nlon_out
781  ! longitude window
782  is = interp%ilon(m,1); ie = interp%ilon(m,2)
783  fis = interp%faci(m,1); fie = interp%faci(m,2)
784  npass = 1
785  dwtsum = 0.
786  wtsum = 0.
787  arsum = 0.
788 
789  ! wrap-around on input grid
790  ! sum using 2 passes (pass 1: end of input grid)
791  if ( ie < is ) then
792  ie = nlon_in
793  fie = 1.0
794  npass = 2
795  endif
796 
797  do np = 1, npass
798  ! pass 2: beginning of input grid
799  if ( np == 2 ) then
800  is = 1
801  fis = 1.0
802  ie = interp%ilon(m,2)
803  fie = interp%faci(m,2)
804  endif
805 
806  ! summing data*weight and weight for single grid point
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) )
813  else
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 )
816  endif
817  enddo
818 
819  if (wtsum > eps) then
820  data_out(m,n) = dwtsum/wtsum
821  if (present(mask_out)) mask_out(m,n) = wtsum/arsum
822  else
823  data_out(m,n) = 0.
824  if (present(mask_out)) mask_out(m,n) = 0.0
825  endif
826 
827  enddo
828  enddo
829 
830  !***********************************************************************
831  ! compute statistics: minimum, maximum, and mean
832  !-----------------------------------------------------------------------
833 
834  if (iverbose > 0) then
835 
836  ! compute statistics of input data
837 
838  call stats(data_in, interp%area_src, asum, dsum, wsum, min_in, max_in, miss_in, mask_in)
839  ! diagnostic messages
840  ! on the root_pe, we can calculate the global mean, minimum and maximum.
841  if(pe == root_pe) then
842  if (wsum > 0.0) then
843  avg_in=dsum/wsum
844  else
845  print *, 'horiz_interp stats: input area equals zero '
846  avg_in=0.0
847  endif
848  if (iverbose > 1) print '(2f16.11)', 'global sum area_in = ', asum, wsum
849  endif
850 
851  ! compute statistics of output data
852  call stats(data_out, interp%area_dst, asum, dsum, wsum, min_out, max_out, miss_out, mask_out)
853  ! diagnostic messages
854  if(pe == root_pe) then
855  if (wsum > 0.0) then
856  avg_out=dsum/wsum
857  else
858  print *, 'horiz_interp stats: output area equals zero '
859  avg_out=0.0
860  endif
861  if (iverbose > 1) print '(2f16.11)', 'global sum area_out = ', asum, wsum
862  endif
863  !---- output statistics ----
864  ! the global mean, min and max are calculated on the root pe.
865  if(pe == root_pe) then
866  write (*,900)
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
871  endif
872 
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)
877 
878  endif
879 
880  !-----------------------------------------------------------------------
881  end subroutine horiz_interp_conserve_version1
882 
883  !#############################################################################
884  subroutine horiz_interp_conserve_version2 ( Interp, data_in, data_out, verbose )
885  !-----------------------------------------------------------------------
886  type(horiz_interp_type), intent(in) :: interp
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
891 
892  data_out = 0.0
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)
897  end do
898 
899  end subroutine horiz_interp_conserve_version2
900 
901  !#######################################################################
902  ! <SUBROUTINE NAME="horiz_interp_conserve_del">
903 
904  ! <OVERVIEW>
905  ! Deallocates memory used by "horiz_interp_type" variables.
906  ! Must be called before reinitializing with horiz_interp_new.
907  ! </OVERVIEW>
908  ! <DESCRIPTION>
909  ! Deallocates memory used by "horiz_interp_type" variables.
910  ! Must be called before reinitializing with horiz_interp_new.
911  ! </DESCRIPTION>
912  ! <TEMPLATE>
913  ! call horiz_interp_conserve_del ( Interp )
914  ! </TEMPLATE>
915 
916  ! <INOUT NAME="Interp" TYPE="horiz_interp_type">
917  ! A derived-type variable returned by previous call
918  ! to horiz_interp_new. The input variable must have
919  ! allocated arrays. The returned variable will contain
920  ! deallocated arrays.
921  ! </INOUT>
922 
923  subroutine horiz_interp_conserve_del ( Interp )
925  type(horiz_interp_type), intent(inout) :: interp
926 
927  select case(interp%version)
928  case (1)
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)
935  case (2)
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)
941  end select
942 
943  end subroutine horiz_interp_conserve_del
944  ! </SUBROUTINE>
945 
946  !#######################################################################
947  !---This statistics is for conservative scheme
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(:,:)
953 
954  integer :: pe, root_pe, npes, p, buffer_int(1)
955  real :: buffer_real(5)
956 
957  pe = mpp_pe()
958  root_pe = mpp_root_pe()
959  npes = mpp_npes()
960 
961  ! sum data, data*area; and find min,max on each pe.
962 
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)
970  else
971  asum = sum(area(:,:))
972  dsum = sum(area(:,:)*dat(:,:))
973  wsum = sum(area(:,:))
974  miss = 0
975  low = minval(dat(:,:))
976  high = maxval(dat(:,:))
977  endif
978 
979  ! other pe send local min, max, avg to the root pe and
980  ! root pe receive these information
981 
982  if(pe == root_pe) then
983  do p = 1, npes - 1
984  ! Force use of "scalar", integer pointer mpp interface
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)
993  enddo
994  else
995  buffer_real(1) = asum
996  buffer_real(2) = dsum
997  buffer_real(3) = wsum
998  buffer_real(4) = low
999  buffer_real(5) = high
1000  ! Force use of "scalar", integer pointer mpp interface
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)
1004  endif
1005 
1006  call mpp_sync_self()
1007 
1008  end subroutine stats
1009 
1010  !#######################################################################
1011 
1012  subroutine data_sum( data, area, facis, facie, facjs, facje, &
1013  dwtsum, wtsum, arsum, mask )
1015  ! sums up the data and weights for a single output grid box
1016  !-----------------------------------------------------------------------
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(:,:)
1021 
1022  ! fac__ = fractional portion of each boundary grid box included
1023  ! in the integral
1024  ! dwtsum = sum(data*area*mask)
1025  ! wtsum = sum(area*mask)
1026  ! arsum = sum(area)
1027  !-----------------------------------------------------------------------
1028  real, dimension(size(area,1),size(area,2)) :: wt
1029  real :: asum
1030  integer :: id, jd
1031  !-----------------------------------------------------------------------
1032 
1033  id=size(area,1); jd=size(area,2)
1034 
1035  wt=area
1036  wt( 1,:)=wt( 1,:)*facis
1037  wt(id,:)=wt(id,:)*facie
1038  wt(:, 1)=wt(:, 1)*facjs
1039  wt(:,jd)=wt(:,jd)*facje
1040 
1041  asum = sum(wt)
1042  arsum = arsum + asum
1043 
1044  if (present(mask)) then
1045  wt = wt * mask
1046  dwtsum = dwtsum + sum(wt*data)
1047  wtsum = wtsum + sum(wt)
1048  else
1049  dwtsum = dwtsum + sum(wt*data)
1050  wtsum = wtsum + asum
1051  endif
1052  !-----------------------------------------------------------------------
1053 
1054  end subroutine data_sum
1055 
1056 
1057  !#######################################################################
1058 
1059 end module horiz_interp_conserve_mod
Definition: fms.F90:20
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()
Definition: fms_io.F90:8566
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)
Definition: mpp.F90:39
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
subroutine horiz_interp_conserve_new_1dx2d(Interp, lon_in, lat_in, lon_out, lat_out, mask_in, mask_out, verbose)
subroutine horiz_interp_conserve_version2(Interp, data_in, data_out, verbose)
#define max(a, b)
Definition: mosaic_util.h:33
subroutine horiz_interp_conserve_version1(Interp, data_in, data_out, verbose, mask_in, mask_out)
#define min(a, b)
Definition: mosaic_util.h:32
real(fp), parameter, public pi