fix: improve modern diag manager performance (#1634)
[FMS.git] / test_fms / diag_manager / test_flush_nc_file.F90
blobb0134f4a0575e08d93a29099f1793d635e063fcf
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 when flush_nc_file=.true.
21 program test_flush_nc_file
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   logical                            :: used             !< Dummy argument to send_data
42   real,                  allocatable :: x(:)             !< X axis data
43   real,                  allocatable :: y(:)             !< Y axis_data
44   real,                  allocatable :: var1_data(:,:)   !< Data for variable 1
45   integer                            :: i                !< For do loops
47   call fms_init
48   call set_calendar_type(JULIAN)
49   call diag_manager_init
51   nx = 360
52   ny = 180
54   allocate(x(nx), y(ny))
55   allocate(var1_data(nx,ny))
56   do i=1,nx
57     x(i) = i
58   enddo
59   do i=1,ny
60     y(i) = -91 + i
61   enddo
63   Time = set_date(2,1,1,0,0,0)
64   Time_step = set_time (3600,0) !< 1 hour
66   id_x  = diag_axis_init('x',  x,  'point_E', 'x', long_name='point_E')
67   id_y  = diag_axis_init('y',  y,  'point_N', 'y', long_name='point_N')
68   id_var1 = register_diag_field  ('atmos', 'ua', (/id_x, id_y/), Time)
70   call diag_manager_set_time_end(set_date(2,1,2,0,0,0))
71   do i = 1, 24
72     Time = Time + Time_step
73     var1_data = real(i)
74     used = send_data(id_var1, var1_data, Time)
75     if (mpp_pe() .eq. mpp_root_pe()) print *, "Calling send_data::", i
76     call diag_send_complete(Time_step)
77     call mpp_sync()
78     ! The file should have been flushed by now
79     call check_answers(i)
80   enddo
82   call diag_manager_end(Time)
83   !call check_answers(i)
84   call fms_end
86   contains
87     subroutine check_answers(time_level)
88       integer, intent(in) :: time_level
89       type(FmsNetcdfFile_t) :: fileobj
90       integer                            :: ntimes
91       real,   allocatable                :: vardata(:,:)
93       if (mpp_pe() .ne. mpp_root_pe()) return
94       if (.not. open_file(fileobj, "test_flush.nc", "read")) &
95         call mpp_error(FATAL, "unable to open test_flush.nc for reading")
97       call get_dimension_size(fileobj, "time", ntimes)
98       if (ntimes .ne. time_level) call mpp_error(FATAL, "time is not the correct size::", time_level)
100       allocate(vardata(nx,ny))
102       call read_data(fileobj, "ua", vardata, unlim_dim_level=i)
103       if (any(vardata .ne. time_level)) &
104         call mpp_error(FATAL, "ua is not the expected result")
106       call close_file(fileobj)
107     end subroutine check_answers
108 end program test_flush_nc_file