FMS  2024.03
Flexible Modeling System
fms_diag_time_reduction.F90
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 !***********************************************************************
19 
20 !> @defgroup fms_diag_time_reduction_mod fms_diag_time_reduction_mod
21 !> @ingroup diag_manager
22 !> @brief fms_diag_time_reduction_mod defines classes encapsulating the diag_manager
23 !! time redution types.
24 !!
25 !> @author Miguel Zuniga
26 !!
27 !> @file
28 !> @brief File for @ref fms_diag_time_reduction_mod
29 !> @addtogroup fms_diag_time_reduction_mod
30 !> @{
31 MODULE fms_diag_time_reduction_mod
32 
33  USE diag_data_mod, only: every_time
34  USE fms_mod, ONLY: error_mesg, fatal
35 
36  implicit none
37 
38  !!TODO: (Future effort) Note that time_diurnal processing is a little different
39  !! and more complex than the other reduction methods, and therefore refactoring its
40  !! processing may simplify the overall related codebase. The refactoring,
41  !! if possible, may be done elsewhere in the diag_manager.
42 
43  !!These parametes are the possible kinds of time reduction operations.
44  INTEGER, PARAMETER :: time_none = 0 !< There is no reduction method
45  INTEGER, PARAMETER :: time_average = 1 !< The reduction method is average
46  INTEGER, PARAMETER :: time_rms = 2 !< The reduction method is rms
47  INTEGER, PARAMETER :: time_max = 3 !< The reduction method is max
48  INTEGER, PARAMETER :: time_min = 4 !< The reduction method is min
49  INTEGER, PARAMETER :: time_sum = 5 !< The reudction method is sum
50  INTEGER, PARAMETER :: time_diurnal = 6 !< The reduction method is diurnal
51  INTEGER, PARAMETER :: time_power = 7 !< The reduction method is power
52 
53  !> @brief Class fmsDiagTimeReduction_type has an encapsulation of the "Fortran enum" time
54  !! reduction integer parameters, plus an encapsulation of the groupings of
55  !! the time reduction types. It is intended to provide some of the functionality
56  !! that was coded in the legacy function diag_data.F90:init_output_fields.
57  !! The functionality in the end is used by send_data in (EFFICIENT) do loops calling
58  !! the weighting or math functions to update buffers.
59  !! The integer parameters above are the legal time reduction types,
60  !! but they are not necessarily mutually exclusive in some contexts.
61  !!
62  !> @addtogroup fms_diag_time_reduction_mod
64  integer , private :: the_time_reduction !< The time reduction type, as an integer defined above.
65  logical , private :: time_averaging !< Set true iff time_average, time_rms, time_power or time_diurnal is true
66  logical , private :: time_ops !< Set true iff time_min, time_max, time_rms or time_average is true.
67  CONTAINS
68  procedure, public :: do_time_averaging => do_time_averaging_imp
69  procedure, public :: has_time_ops => has_time_ops_imp
70  procedure, public :: is_time_none => is_time_none_imp
71  procedure, public :: is_time_average => is_time_average_imp
72  procedure, public :: is_time_rms => is_time_rms_imp
73  procedure, public :: is_time_max => is_time_max_imp
74  procedure, public :: is_time_min => is_time_min_imp
75  procedure, public :: is_time_sum => is_time_sum_imp
76  procedure, public :: is_time_diurnal => is_time_diurnal_imp
77  procedure, public :: is_time_power => is_time_power_imp
78  procedure, public :: initialize
79  procedure, public :: copy
81 
82  !> @brief This interface is for the class constructor.
83  !> @addtogroup fms_diag_time_reduction_mod
86  end interface fmsdiagtimereduction_type
87 
88 CONTAINS
89 
90  !> @brief The class contructors. Just allocates the class and calls an initializer
91  !! @return An allocated instance of fmsDiagTimeReduction_type, which is nitialized using
92  !! provided values for arguments dt and out_freqeuncy.
93  function fmsdiagtimereduction_type_constructor(dt, out_frequency) result(time_redux)
94  integer, intent(in) :: dt !< The redution type (time_rms, time_power, etc)
95  integer, intent(in) :: out_frequency !< The output frequency.
96  class(fmsdiagtimereduction_type), allocatable :: time_redux !< The instance of the fmsDiagTimeReduction_type
97  !!class allocated and returned by this constructor.
98  allocate(time_redux)
99  call time_redux%initialize(dt, out_frequency)
101 
102  !> @brief Initialize the object. As an alternative to the constructor, one can
103  !! allocate an fmsDiagTimeReduction_type instance, then call its initialize function.
104  subroutine initialize(this, dt, out_frequency)
105  class(fmsdiagtimereduction_type), intent(inout) :: this !< The fmsDiagTimeReduction_type object
106  integer, intent(in) :: dt !< The redution type (time_rms, time_porer, etc)
107  integer, intent(in) :: out_frequency !< The output frequency.
108 
109  this%the_time_reduction = dt
110 
111  !! Set the time_averaging flag
112  !! See legacy init_ouput_fields function, lines 1470ff
113  IF(( dt .EQ. time_average) .OR. (dt .EQ. time_rms) .OR. (dt .EQ. time_power) .OR. &
114  & (dt .EQ. time_diurnal)) THEN
115  this%time_averaging = .true.
116  ELSE
117  this%time_averaging= .false.
118  IF((dt .NE. time_max) .AND. (dt .ne. time_min) .AND. (dt .NE. time_sum) &
119  & .AND. (dt .NE. time_none)) THEN
120  CALL error_mesg('fmsDiagTimeReduction_type: initialize', &
121  & 'time_averaging=.false. but reduction type not compatible', fatal)
122  ENDIF
123  END IF
124 
125  !!TODO: (MDM) Add other checks? E.g. If time_averaging == .false., then
126  !! out_frequency == EVERY_TIME
127 
128  IF((dt .EQ. time_min) .OR. (dt .EQ. time_max) .OR. &
129  & ( dt .EQ. time_average) .OR. (dt .EQ. time_sum) ) THEN
130  this%time_ops = .true.
131  ELSE
132  this%time_ops = .false.
133  END IF
134  end subroutine initialize
135 
136  !> @brief Copy the source time reduction object into the this object.
137  subroutine copy(this, source)
138  class(fmsdiagtimereduction_type),intent(inout):: this !< The fmsDiagTimeReduction_type object
139  class(fmsdiagtimereduction_type),intent(in):: source !< The fmsDiagTimeReduction_type object
140  this%the_time_reduction = source%the_time_reduction
141  this%time_averaging = source%time_averaging
142  this%time_ops = source%time_ops
143  end subroutine copy
144 
145  !> \brief Returns true if any of time_min, time_max, time_rms or time_average is true.
146  !! @return true if any of time_min, time_max, time_rms or time_average is true.
147  pure function has_time_ops_imp (this)
148  class(fmsdiagtimereduction_type), intent(in) :: this !<The object this function is bound to.
149  logical :: has_time_ops_imp
150  has_time_ops_imp = this%time_ops
151  end function has_time_ops_imp
152 
153  !> \brief Returns true iff time_averaging is true.
154  !! @return true iff time_averaging is true.
155  pure function do_time_averaging_imp (this)
156  class(fmsdiagtimereduction_type), intent(in) :: this !<The object this function is bound to.
157  logical :: do_time_averaging_imp
158  do_time_averaging_imp = this%time_averaging
159  end function do_time_averaging_imp
160 
161  !> \brief Returns true iff the_time_reduction is time_average
162  !! @return true iff the_time_reduction is time_average
163  pure function is_time_average_imp (this)
164  class(fmsdiagtimereduction_type), intent(in) :: this !<The object this function is bound to.
165  logical :: is_time_average_imp
166  is_time_average_imp = this%the_time_reduction .EQ. time_average
167  end function is_time_average_imp
168 
169  !> \brief Returns true iff the_time_reduction is time_none
170  !! @return true iff the_time_reduction is time_none
171  pure function is_time_none_imp (this)
172  class(fmsdiagtimereduction_type), intent(in) :: this !<The object this function is bound to.
173  logical :: is_time_none_imp
174  is_time_none_imp = (this%the_time_reduction .EQ. time_none)
175  end function is_time_none_imp
176 
177  !> \brief Returns true iff the_time_reduction is time_rms
178  !! @return true iff the_time_reduction is time_rms
179  pure function is_time_rms_imp (this)
180  class(fmsdiagtimereduction_type), intent(in) :: this !<The object this function is bound to.
181  logical :: is_time_rms_imp
182  is_time_rms_imp = this%the_time_reduction .EQ. time_rms
183  end function is_time_rms_imp
184 
185  !> \brief Returns true iff the_time_reduction is time_max
186  !! @return true iff the_time_reduction is time_max
187  pure function is_time_max_imp (this)
188  class(fmsdiagtimereduction_type), intent(in) :: this !<The object this function is bound to.
189  logical :: is_time_max_imp
190  is_time_max_imp = this%the_time_reduction .EQ. time_max
191  end function is_time_max_imp
192 
193  !> \brief Returns true iff the_time_reduction is time_min
194  !! @return true iff the_time_reduction is time_min
195  pure function is_time_min_imp (this)
196  class(fmsdiagtimereduction_type), intent(in) :: this !<The object this function is bound to.
197  logical :: is_time_min_imp
198  is_time_min_imp = this%the_time_reduction .EQ. time_min
199  end function is_time_min_imp
200 
201  !> \brief Returns true iff the_time_reduction is time_sum
202  !! @return true iff the_time_reduction is time_sum
203  pure function is_time_sum_imp (this)
204  class(fmsdiagtimereduction_type), intent(in) :: this !<The object this function is bound to.
205  logical :: is_time_sum_imp
206  is_time_sum_imp = this%the_time_reduction .EQ. time_sum
207  end function is_time_sum_imp
208 
209  !> \brief Returns true iff the_time_reduction is time_diurnal
210  !! @return true iff the_time_reduction is time_diurnal
211  pure function is_time_diurnal_imp (this)
212  class(fmsdiagtimereduction_type), intent(in) :: this !<The object this function is bound to.
213  logical :: is_time_diurnal_imp
214  is_time_diurnal_imp = this%the_time_reduction .EQ. time_diurnal
215  end function is_time_diurnal_imp
216 
217  !> \brief Returns true iff the_time_reduction is time_power
218  !! @return true iff the_time_reduction is time_power
219  pure function is_time_power_imp (this)
220  class(fmsdiagtimereduction_type), intent(in) :: this !<The object this function is bound to.
221  logical :: is_time_power_imp
222  is_time_power_imp = this%the_time_reduction .EQ. time_power
223  end function is_time_power_imp
224 
225 end module fms_diag_time_reduction_mod
226 !> @}
227 ! close documentation grouping
pure logical function is_time_none_imp(this)
Returns true iff the_time_reduction is time_none.
pure logical function is_time_sum_imp(this)
Returns true iff the_time_reduction is time_sum.
pure logical function do_time_averaging_imp(this)
Returns true iff time_averaging is true.
pure logical function is_time_min_imp(this)
Returns true iff the_time_reduction is time_min.
class(fmsdiagtimereduction_type) function, allocatable fmsdiagtimereduction_type_constructor(dt, out_frequency)
The class contructors. Just allocates the class and calls an initializer.
pure logical function is_time_rms_imp(this)
Returns true iff the_time_reduction is time_rms.
pure logical function is_time_max_imp(this)
Returns true iff the_time_reduction is time_max.
subroutine copy(this, source)
Copy the source time reduction object into the this object.
integer, parameter time_min
The reduction method is min.
integer, parameter time_diurnal
The reduction method is diurnal.
integer, parameter time_power
The reduction method is power.
pure logical function is_time_power_imp(this)
Returns true iff the_time_reduction is time_power.
pure logical function is_time_diurnal_imp(this)
Returns true iff the_time_reduction is time_diurnal.
integer, parameter time_average
The reduction method is average.
integer, parameter time_sum
The reudction method is sum.
integer, parameter time_rms
The reduction method is rms.
subroutine initialize(this, dt, out_frequency)
Initialize the object. As an alternative to the constructor, one can allocate an fmsDiagTimeReduction...
pure logical function is_time_average_imp(this)
Returns true iff the_time_reduction is time_average.
integer, parameter time_max
The reduction method is max.
pure logical function has_time_ops_imp(this)
Returns true if any of time_min, time_max, time_rms or time_average is true.
subroutine, public error_mesg(routine, message, level)
Print notes, warnings and error messages; terminates program for warning and error messages....
Definition: fms.F90:498