FMS  2024.03
Flexible Modeling System
fms_diag_reduction_methods.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_reduction_methods_mod fms_diag_reduction_methods_mod
21 !> @ingroup diag_manager
22 !! @brief fms_diag_reduction_methods_mod contains routines that are meant to be used for
23 !! error checking and setting up to do the reduction methods
24 
25 !> @file
26 !> @brief File for @ref fms_diag_reduction_methods_mod
27 
28 !> @addtogroup fms_diag_reduction_methods_mod
29 !> @{
30 module fms_diag_reduction_methods_mod
31  use platform_mod, only: r8_kind, r4_kind
32  use fms_diag_bbox_mod, only: fmsdiagibounds_type
33  use fms_string_utils_mod, only: string
34  use diag_data_mod, only: time_diurnal, time_rms
35  use mpp_mod
36  implicit none
37  private
38 
41 
42  !> @brief Does the time_none reduction method. See include/fms_diag_reduction_methods.inc
43  !TODO This needs to be extended to integers
44  interface do_time_none
45  module procedure do_time_none_r4, do_time_none_r8
46  end interface do_time_none
47 
48  !> @brief Does the time_min reduction method. See include/fms_diag_reduction_methods.inc
49  !TODO This needs to be extended to integers
50  interface do_time_min
51  module procedure do_time_min_r4, do_time_min_r8
52  end interface do_time_min
53 
54  !> @brief Does the time_max reduction method. See include/fms_diag_reduction_methods.inc
55  !TODO This needs to be extended to integers
56  interface do_time_max
57  module procedure do_time_max_r4, do_time_max_r8
58  end interface do_time_max
59 
60  !> @brief Sum update updates the buffer for any reductions that involve summation
61  !! (ie. time_sum, avg, rms, pow)
62  !!TODO This needs to be extended to integers
64  module procedure do_time_sum_update_r4, do_time_sum_update_r8
65  end interface
66 
67  !> @brief Finishes a reduction that involves an average
68  !! (ie. time_avg, rms, pow)
69  !! This takes the average at the end of the time step
70  interface time_update_done
71  module procedure sum_update_done_r4, sum_update_done_r8
72  end interface
73 
74  contains
75 
76  !> @brief Checks improper combinations of is, ie, js, and je.
77  !! @return The error message, empty string if no errors were found
78  !> @note accept_data works in either one or another of two modes.
79  !! 1. Input field is a window (e.g. FMS physics)
80  !! 2. Input field includes halo data
81  !! It cannot handle a window of data that has halos.
82  !! (A field with no windows or halos can be thought of as a special case of either mode.)
83  !! The logic for indexing is quite different for these two modes, but is not clearly separated.
84  !! If both the beggining and ending indices are present, then field is assumed to have halos.
85  !! If only beggining indices are present, then field is assumed to be a window.
86  !> @par
87  !! There are a number of ways a user could mess up this logic, depending on the combination
88  !! of presence/absence of is,ie,js,je. The checks below should catch improper combinations.
89  pure function check_indices_order(is_in, ie_in, js_in, je_in) &
90  result(error_msg)
91  integer, intent(in), optional :: is_in, ie_in, js_in, je_in !< Indices passed to fms_diag_accept_data()
92  character(len=128) :: error_msg !< An error message used only for testing purpose!!!
93 
94  error_msg = ""
95  IF ( PRESENT(ie_in) ) THEN
96  IF ( .NOT.PRESENT(is_in) ) THEN
97  error_msg = 'ie_in present without is_in'
98  return
99  END IF
100  IF ( PRESENT(js_in) .AND. .NOT.PRESENT(je_in) ) THEN
101  error_msg = 'is_in and ie_in present, but js_in present without je_in'
102  return
103  END IF
104  END IF
105 
106  IF ( PRESENT(je_in) ) THEN
107  IF ( .NOT.PRESENT(js_in) ) THEN
108  error_msg = 'je_in present without js_in'
109  return
110  END IF
111  IF ( PRESENT(is_in) .AND. .NOT.PRESENT(ie_in) ) THEN
112  error_msg = 'js_in and je_in present, but is_in present without ie_in'
113  return
114  END IF
115  END IF
116  end function check_indices_order
117 
118  !> @brief Sets the logical mask based on mask or rmask
119  !> @return logical mask
120  function init_mask(rmask, mask, field) &
121  result(oor_mask)
122  LOGICAL, DIMENSION(:,:,:,:), allocatable, INTENT(in) :: mask !< The location of the mask
123  CLASS(*), DIMENSION(:,:,:,:), allocatable, INTENT(in) :: rmask !< The masking values
124  CLASS(*), DIMENSION(:,:,:,:), intent(in) :: field !< Field_data
125 
126  logical, allocatable, dimension(:,:,:,:) :: oor_mask !< mask
127 
128  ALLOCATE(oor_mask(SIZE(field, 1), SIZE(field, 2), SIZE(field, 3), SIZE(field, 4)))
129  oor_mask = .true.
130 
131  if (allocated(mask)) then
132  oor_mask = mask
133  elseif (allocated(rmask)) then
134  select type (rmask)
135  type is (real(kind=r8_kind))
136  WHERE (rmask < 0.5_r8_kind) oor_mask = .false.
137  type is (real(kind=r4_kind))
138  WHERE (rmask < 0.5_r4_kind) oor_mask = .false.
139  end select
140  endif
141 
142  end function init_mask
143 
144  !> @brief Sets the weight based on the weight passed into send_data (1.0_r8_kind if the weight is not passed in)
145  !! The weight will be saved as an r8 and converted to r4 as needed
146  !! @return weight to use when averaging
147  pure function set_weight(weight) &
148  result(out_weight)
149  CLASS(*), INTENT(in), OPTIONAL :: weight !< The weight use when averaging
150 
151  real(kind=r8_kind) :: out_weight
152 
153  out_weight = 1.0_r8_kind
154  if (present(weight)) then
155  select type(weight)
156  type is (real(kind=r8_kind))
157  out_weight = real(weight, kind = r8_kind)
158  type is (real(kind=r4_kind))
159  out_weight = real(weight, kind = r8_kind)
160  end select
161  endif
162  end function set_weight
163 
164 #include "fms_diag_reduction_methods_r4.fh"
165 #include "fms_diag_reduction_methods_r8.fh"
166 
167 end module fms_diag_reduction_methods_mod
168 !> @}
169 ! close documentation grouping
integer, parameter time_diurnal
The reduction method is diurnal.
Definition: diag_data.F90:122
integer, parameter time_rms
The reudction method is root mean square of values.
Definition: diag_data.F90:121
pure real(kind=r8_kind) function, public set_weight(weight)
Sets the weight based on the weight passed into send_data (1.0_r8_kind if the weight is not passed in...
logical function, dimension(:,:,:,:), allocatable, public init_mask(rmask, mask, field)
Sets the logical mask based on mask or rmask.
pure character(len=128) function, public check_indices_order(is_in, ie_in, js_in, je_in)
Checks improper combinations of is, ie, js, and je.
Does the time_max reduction method. See include/fms_diag_reduction_methods.inc.
Does the time_min reduction method. See include/fms_diag_reduction_methods.inc.
Sum update updates the buffer for any reductions that involve summation (ie. time_sum,...
Finishes a reduction that involves an average (ie. time_avg, rms, pow) This takes the average at the ...
character(:) function, allocatable, public string(v, fmt)
Converts a number or a Boolean value to a string.
Does the time_none reduction method. See include/fms_diag_reduction_methods.inc.
Data structure holding a 3D bounding box. It is commonlyused to represent the interval bounds or limi...