chore: append -dev to version number (#1641)
[FMS.git] / test_fms / diag_manager / test_prepend_date.F90
blob24a5ae2986a2877114d6c01a0699fb2183395dcc
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 init date is prepended to the file name
21 program test_prepend_date
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, input_nml_file
29   use fms2_io_mod,      only: FmsNetcdfFile_t, open_file, close_file, read_data, get_dimension_size
30   use platform_mod,     only: r4_kind
32   implicit none
34   integer         :: id_var0, id_var2, id_var1 !< diag field ids
35   integer         :: id_axis1         !< Id for axis
36   logical         :: used             !< for send_data calls
37   integer         :: ntimes = 48      !< Number of time steps
38   real            :: vdata            !< Buffer to store the data
39   type(time_type) :: Time             !< "Model" time
40   type(time_type) :: Time_step        !< Time step for the "simulation"
41   integer         :: i                !< For do loops
42   logical         :: pass_diag_time = .True.   !< .True. if passing the time to diag_manager_init
44   integer :: io_status !< Status when reading the namelist
46   namelist / test_prepend_date_nml / pass_diag_time
48   call fms_init
50   read (input_nml_file, test_prepend_date_nml, iostat=io_status)
51   if (io_status > 0) call mpp_error(FATAL,'=>test_prepend_date: Error reading input.nml')
53   call set_calendar_type(JULIAN)
55   ! This is going to be different from the base_date
56   if (pass_diag_time) then
57     call diag_manager_init(time_init=(/2, 1, 1, 0, 0, 0/))
58   else
59     call diag_manager_init()
60   endif
62   Time = set_date(2,1,1,0,0,0)
63   Time_step = set_time (3600,0) !< 1 hour
64   call diag_manager_set_time_end(set_date(2,1,3,0,0,0))
66   id_axis1 = diag_axis_init('dummy_axis', (/real(1.)/), "mullions", "X")
67   id_var0 = register_diag_field  ('ocn_mod', 'var0', Time)
68   id_var2 = register_static_field ('ocn_mod', 'var2', (/id_axis1/))
70   ! This is a different start_time, should lead to a crash if the variable is in the diag table yaml
71   id_var1 = register_diag_field  ('ocn_mod', 'var1', set_date(2,1,6,0,0,0))
73   used = send_data(id_var2, real(123.456))
74   do i = 1, ntimes
75     Time = Time + Time_step
76     vdata = real(i)
78     used = send_data(id_var0, vdata, Time) !< Sending data every hour!
80     call diag_send_complete(Time_step)
81   enddo
83   call diag_manager_end(Time)
85   call check_output()
86   call fms_end
88   contains
90   !< @brief Check the diag manager output
91   subroutine check_output()
92     type(FmsNetcdfFile_t) :: fileobj     !< Fms2io fileobj
93     integer               :: var_size    !< Size of the variable reading
94     real(kind=r4_kind), allocatable     :: var_data(:) !< Buffer to read variable data to
95     integer               :: j           !< For looping
97     if (.not. open_file(fileobj, "00020101.test_non_static.nc", "read")) &
98       call mpp_error(FATAL, "Error opening file:00020101.test_non_static.nc to read")
100     call get_dimension_size(fileobj, "time", var_size)
101     if (var_size .ne. 48) call mpp_error(FATAL, "The dimension of time in the file:test_0days is not the "//&
102                                                  "correct size!")
103     allocate(var_data(var_size))
104     var_data = -999.99
106     call read_data(fileobj, "var0", var_data)
107     do j = 1, var_size
108       if (var_data(j) .ne. real(j, kind=r4_kind)) call mpp_error(FATAL, "The variable data for var1 at time level:"//&
109                                                           string(j)//" is not the correct value!")
110     enddo
112     call close_file(fileobj)
114     if (.not. open_file(fileobj, "00020101.test_static.nc", "read")) &
115       call mpp_error(FATAL, "Error opening file:00020101.test_static.nc to read")
117     call read_data(fileobj, "var2", var_data(1))
118     if (var_data(1) .ne. real(123.456, kind=r4_kind)) call mpp_error(FATAL, &
119       "The variable data for var2 is not the correct value!")
121     call close_file(fileobj)
123   end subroutine check_output
124 end program test_prepend_date