chore: append -dev to version number (#1641)
[FMS.git] / test_fms / diag_manager / test_output_every_freq.F90
blob3cc88d08acf4696c0d33c2121b7a6e8ea9a09cce
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 diag manager when the file frequency is set to 0 days
21 program test_output_every_freq
23   use fms_mod,          only: fms_init, fms_end, string
24   use diag_manager_mod, only: diag_axis_init, send_data, diag_send_complete, diag_manager_set_time_end, &
25                               register_diag_field, diag_manager_init, diag_manager_end, register_static_field, &
26                               diag_axis_init
27   use time_manager_mod, only: time_type, operator(+), JULIAN, set_time, set_calendar_type, set_date
28   use mpp_mod,          only: FATAL, mpp_error
29   use fms2_io_mod,      only: FmsNetcdfFile_t, open_file, close_file, read_data, get_dimension_size
31   implicit none
33   integer         :: id_var0, id_var1, id_var2 !< diag field ids
34   integer         :: id_axis1         !< Id for axis
35   logical         :: used             !< for send_data calls
36   integer         :: ntimes = 48      !< Number of time steps
37   real            :: vdata            !< Buffer to store the data
38   type(time_type) :: Time             !< "Model" time
39   type(time_type) :: Time_step        !< Time step for the "simulation"
40   integer         :: i                !< For do loops
42   call fms_init
43   call set_calendar_type(JULIAN)
44   call diag_manager_init
46   Time = set_date(2,1,1,0,0,0)
47   Time_step = set_time (3600,0) !< 1 hour
48   call diag_manager_set_time_end(set_date(2,1,3,0,0,0))
50   id_axis1 = diag_axis_init('dummy_axis', (/real(1.)/), "mullions", "X")
51   id_var0 = register_diag_field  ('ocn_mod', 'var0', Time)
52   id_var1 = register_diag_field  ('ocn_mod', 'var1', Time)
53   id_var2 = register_static_field ('ocn_mod', 'var2', (/id_axis1/))
55   used = send_data(id_var2, real(123.456))
56   do i = 1, ntimes
57     Time = Time + Time_step
58     vdata = real(i)
60     used = send_data(id_var0, vdata, Time) !< Sending data every hour!
61     if (mod(i,2) .eq. 0) used = send_data(id_var1, vdata, Time) !< Sending data every two hours!
63     call diag_send_complete(Time_step)
64   enddo
66   call diag_manager_end(Time)
68   call check_output()
69   call fms_end
71   contains
73   !< @brief Check the diag manager output
74   subroutine check_output()
75     type(FmsNetcdfFile_t) :: fileobj     !< Fms2io fileobj
76     integer               :: var_size    !< Size of the variable reading
77     real, allocatable     :: var_data(:) !< Buffer to read variable data to
78     integer               :: j           !< For looping
80     if (.not. open_file(fileobj, "test_0days.nc", "read")) &
81       call mpp_error(FATAL, "Error opening file:test_0days.nc to read")
83     call get_dimension_size(fileobj, "time", var_size)
84     if (var_size .ne. 48) call mpp_error(FATAL, "The dimension of time in the file:test_0days is not the "//&
85                                                  "correct size!")
86     allocate(var_data(var_size))
87     var_data = -999.99
89     call read_data(fileobj, "var0", var_data)
90     do j = 1, var_size
91       if (var_data(j) .ne. real(j)) call mpp_error(FATAL, "The variable data for var1 at time level:"//&
92                                                           string(j)//" is not the correct value!")
93     enddo
95     call close_file(fileobj)
96   end subroutine check_output
97 end program test_output_every_freq