fix: improve modern diag manager performance (#1634)
[FMS.git] / test_fms / diag_manager / test_multiple_send_data.F90
blobd8df31b42f2e6611f780987bcb1faa98a3dcdc22
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.
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.
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 !> @brief  This programs tests fields that call send_data multiple times
21 program test_multiple_send_data
22   use fms_mod, only: fms_init, fms_end
23   use diag_manager_mod
24   use mpp_mod
25   use mpp_domains_mod
26   use platform_mod,     only: r8_kind, r4_kind
27   use time_manager_mod, only: time_type, set_calendar_type, set_date, JULIAN, set_time, OPERATOR(+)
28   use fms2_io_mod
29   use fms_diag_yaml_mod
31   implicit none
33   type(time_type)                    :: Time             !< Time of the simulation
34   type(time_type)                    :: Time_step        !< Time_step of the simulation
35   integer                            :: nx               !< Number of x points
36   integer                            :: ny               !< Number of y points
37   integer                            :: nz               !< Number of z points
38   integer                            :: id_x             !< Axis id for the x dimension
39   integer                            :: id_y             !< Axis id for the y dimension
40   integer                            :: id_var1          !< Field id for 1st variable
41   integer                            :: id_var2          !< Field id for 2nd variable
42   logical                            :: used             !< Dummy argument to send_data
43   real,                  allocatable :: x(:)             !< X axis data
44   real,                  allocatable :: y(:)             !< Y axis_data
45   real,                  allocatable :: var1_data(:,:)   !< Data for variable 1
46   logical,               allocatable :: var1_mask(:,:)   !< Mask for variable 1
47   integer                            :: i                !< For do loops
49   call fms_init
50   call set_calendar_type(JULIAN)
51   call diag_manager_init
53   nx = 360
54   ny = 180
56   allocate(x(nx), y(ny))
57   allocate(var1_data(nx,ny), var1_mask(nx,ny))
58   do i=1,nx
59     x(i) = i
60   enddo
61   do i=1,ny
62     y(i) = -91 + i
63   enddo
65   Time = set_date(2,1,1,0,0,0)
66   Time_step = set_time (3600,0) !< 1 hour
68   id_x  = diag_axis_init('x',  x,  'point_E', 'x', long_name='point_E')
69   id_y  = diag_axis_init('y',  y,  'point_N', 'y', long_name='point_N')
71   ! id_var1 is using diag manager similarly to how the `rv_ice` and `rv_T` in the river code
72   id_var1 = register_diag_field  ('atmos', 'ua', (/id_x, id_y/), Time, missing_value=-999., mask_variant=.True., &
73     multiple_send_data=.True.)
75   ! id_var2 is using diag manager similarly to way it is used in the vert_diff module code
76   id_var2 = register_diag_field  ('atmos', 'va', (/id_x, id_y/), Time, missing_value=-999., &
77     multiple_send_data=.True.)
79   call diag_manager_set_time_end(set_date(2,1,2,0,0,0))
80   do i = 1, 24
81     Time = Time + Time_step
83     var1_data = real(i)
85     var1_mask = .False.
86     ! Only count the data for the (:,1) section of the grid on this send_data call
87     var1_mask(:,1) = .True.
88     used = send_data(id_var1, var1_data, Time, mask=var1_mask)
90     var1_mask = .False.
91     ! Only count the data for the (:,2:) section of the grid on this send_data call
92     var1_mask(:,2:) = .True.
93     used = send_data(id_var1, var1_data, Time, mask=var1_mask)
95     used = send_data(id_var2, var1_data, Time)
96     used = send_data(id_var2, var1_data, Time)
97     call diag_send_complete(Time_step)
98   enddo
100   call diag_manager_end(Time)
101   call check_answers()
102   call fms_end
104   contains
105     subroutine check_answers()
106       type(FmsNetcdfFile_t) :: fileobj
107       integer                            :: ntimes
108       integer                            :: nx
109       integer                            :: ny
110       real,   allocatable                :: vardata(:,:)
111       real                               :: ans_var
112       integer                            :: i, j
114       if (.not. open_file(fileobj, "test_multiple_sends.nc", "read")) &
115         call mpp_error(FATAL, "unable to open test_var_masks.nc for reading")
117       call get_dimension_size(fileobj, "time", ntimes)
118       if (ntimes .ne. 1) call mpp_error(FATAL, "time is not the correct size!")
120       call get_dimension_size(fileobj, "x", nx)
121       if (nx .ne. 360) call mpp_error(FATAL, "x is not the correct size!")
123       call get_dimension_size(fileobj, "y", ny)
124       if (ny .ne. 180) call mpp_error(FATAL, "y is not the correct size!")
126       allocate(vardata(nx,ny))
128       ans_var = 0.
129       do i = 1, 24
130         ans_var = ans_var + real(i)
131       enddo
132       ans_var = ans_var / 24
134       call read_data(fileobj, "ua", vardata)
135       if (any(vardata .ne. ans_var)) &
136         call mpp_error(FATAL, "ua is not the expected result")
138       call read_data(fileobj, "va", vardata)
139       if (any(vardata .ne. ans_var)) &
140         call mpp_error(FATAL, "va is not the expected result")
142       call close_file(fileobj)
143     end subroutine check_answers
144 end program test_multiple_send_data