FMS  2025.04
Flexible Modeling System
fms_diag_outfield.F90
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 
19 !> @defgroup fms_diag_outfield_mod fms_diag_outfield_mod
20 !> @ingroup diag_manager
21 !> @brief fms_diag_outfield_mod defines data types and utility or auxiliary routines
22 !! useful in updating the output buffer.
23 !!
24 !> @author Miguel Zuniga
25 !!
26 !! <TT>fms_diag_outfield_mod</TT> The output buffer updating routines are passed configuration
27 !! and control data with types defined in this module; and some utility functions called by the
28 !! updating routines are
29 !! defined here.
30 !!
31 !> @file
32 !> @brief File for @ref fms_diag_outfield_mod
33 !> @addtogroup fms_diag_outfield_mod
34 !> @{
35 MODULE fms_diag_outfield_mod
36  USE platform_mod
37  USE mpp_mod, only :fatal, warning
38  USE fms_mod, only :lowercase, uppercase, error_mesg, fms_error_handler
39 
40 
41  !! TODO: these might need removal or replacement
42  USE diag_data_mod, only:time_zero
43  USE diag_data_mod, only: glo_reg_val, glo_reg_val_alt, region_out_use_alt_value, very_large_axis_length, coord_type
44  USE diag_data_mod, only: fmsdiagibounds_type, input_field_type, output_field_type
45  USE fms_diag_time_reduction_mod, only: fmsdiagtimereduction_type, time_none , time_average, time_rms
46  USE fms_diag_time_reduction_mod, only: time_max, time_min, time_sum, time_power
47 
48  implicit none
49 
50  !> @brief Class fmsDiagOutfield_type (along with class ms_diag_outfield_index_type )
51  !! contain information used in updating the output buffers by the diag_manager
52  !! send_data routines. In some sense they can be seen as encapsulating related
53  !! information in a convenient way (e.g. to pass to functions and for do loop
54  !! controls.)
55  !!
56  !! Class fmsDiagOutfield_type also contains a significant subset of the fields
57  !! and routines of of the legacy class output_field_type
58  !! TODO: (MDM) This class will need further development for the modern_diag effort.
59  !! For its development, consider the legacy diag_util::init_output_field already
60  !! in place. Fields added so are used the field buffer math/dmUpdate functions.
61  !! TODO (MDM) : Should the MDM have pow_value be type REAL?
62  !> @ingroup fms_diag_outfield_mod
63  TYPE, public :: fmsdiagoutfield_type
64  PRIVATE
65  CHARACTER(len=:), ALLOCATABLE :: module_name !< Module name.
66  CHARACTER(len=:), ALLOCATABLE :: field_name !< Output field name.
67  CHARACTER(len=:), ALLOCATABLE :: output_name !< Output name written to file.
68  CHARACTER(len=:), ALLOCATABLE :: output_file !< File where field should be written.
69 
70  !!Major outer loop controls in send_data functions.
71  INTEGER :: pow_value !< Power value for rms or pow(x) calculations
72  LOGICAL :: phys_window !< TODO: Rename? OMP subsetted data, See output_fields
73  LOGICAL :: need_compute !< True iff is local_output and current PE take part in send_data.
74  LOGICAL :: reduced_k_range !< If true, the local start and end indecies are used in k (i.e. 3rd) dim.
75  LOGICAL :: missvalue_present !<
76  LOGICAL :: mask_variant
77  LOGICAL :: mask_present !< True iff mask argument is present in user-facing send function call.
78  !< Note this field exists since the actual mask argument in the send
79  !< function call may be downstream replaced by a null pointer which
80  !< is considered present.
81 
82  TYPE(fmsdiagtimereduction_type) :: time_reduction !< Instance of the fmsDiagTimeTeduction_type.
83 
84  !!TODO (Future effort? ) : a pointer for time_min and time_max comparison function
85  !! If possible, this can remove the innermost if/then/else construct in the buffer update loops.
86  !! min_max_f_ptr => (should point to < or > operators)
87 
88  !! gcc error: Interface ‘addwf’ at (1) must be explicit
89  ! procedure (addwf), pointer, nopass :: f_ptr => null () !!A pointer to the field weighing procedure
90 
91  CONTAINS
92  procedure :: get_module_name
93  procedure :: get_field_name
94  procedure :: get_output_name
95  procedure :: get_output_file
96  procedure :: get_pow_value
97  procedure :: get_phys_window
98  procedure :: get_need_compute
99  procedure :: get_reduced_k_range
100  procedure :: get_missvalue_present
101  procedure :: get_mask_variant
102  procedure :: get_mask_present
103  procedure :: get_time_reduction
104  procedure, public :: initialize => initialize_outfield_imp
105  procedure :: initialize_for_ut
106 
107  END TYPE fmsdiagoutfield_type
108 
109  !> @brief Class fms_diag_outfield_index_type which (along with class fmsDiagOutfield_type)
110  !! encapsulate related information used in updating the output buffers by the diag_manager
111  !! send_data routines. This class in particular focuses on do loop index controls or settings.
112  !! Note that the index names in this class should be indentical to the names used in the
113  !! diag_manager send_data functions and in the "math" buffer update functions. The purpose
114  !! of this class is also to allow for a smaller call function signature for the math/buffer
115  !! update functions.
116  !> @ingroup fms_diag_outfield_mod
118  PRIVATE
119  INTEGER :: f1,f2 !< Indecies used specify 1st dim bounds of field, mask and rmask.
120  INTEGER :: f3,f4 !< Indecies used specify 2st dim bounds of field, mask and rmask.
121  INTEGER :: is, js, ks !< Start indecies in each spatial dim of the field_data; and
122  !! may be user provided in send_data
123  Integer :: ie, je, ke !< End indecies in each spatial dim of the field_data; and
124  !! may be user provided in send_data
125  INTEGER :: hi !< halo size in x direction. Same name as in send_data
126  INTEGER :: hj !< halo size in y direction. Same
127  CONTAINS
128  procedure :: initialize => initialize_outfield_index_type
129  procedure :: get_f1
130  procedure :: get_f2
131  procedure :: get_f3
132  procedure :: get_f4
133  procedure :: get_is
134  procedure :: get_js
135  procedure :: get_ks
136  procedure :: get_ie
137  procedure :: get_je
138  procedure :: get_ke
139  procedure :: get_hi
140  procedure :: get_hj
142 
143 CONTAINS
144 
145  !> @brief Gets module_name
146  !! @return copy of the module_name character array
147  pure function get_module_name (this) &
148  result(rslt)
149  class(fmsdiagoutfield_type), intent(in) :: this !< diag object
150  character(len=:), allocatable :: rslt
151  rslt = this%module_name
152  end function get_module_name
153 
154  !> @brief Gets field_name
155  !! @return copy of the field_name character array
156  pure function get_field_name (this) &
157  result(rslt)
158  class(fmsdiagoutfield_type), intent(in) :: this !< diag object
159  character(len=:), allocatable :: rslt
160  rslt = this%field_name
161  end function get_field_name
162 
163  !> @brief Gets output_name
164  !! @return copy of the output_name character array
165  pure function get_output_name (this) &
166  result(rslt)
167  class(fmsdiagoutfield_type), intent(in) :: this !< diag object
168  character(len=:), allocatable :: rslt
169  rslt = this%output_name
170  end function get_output_name
171 
172  !> @brief Gets output_file
173  !! @return copy of the output_file character array
174  pure function get_output_file (this) &
175  result(rslt)
176  class(fmsdiagoutfield_type), intent(in) :: this !< diag object
177  character(len=:), allocatable :: rslt
178  rslt = this%output_file
179  end function get_output_file
180 
181  !> @brief Gets pow_value
182  !! @return copy of integer member pow_value
183  pure integer function get_pow_value (this) result(rslt)
184  class(fmsdiagoutfield_type), intent(in) :: this !< The fmsDiagOutfield_type object
185  rslt = this%pow_value
186  end function get_pow_value
187 
188  !> @brief Gets phys_window
189  !! @return copy of integer member phys_window
190  pure logical function get_phys_window (this) result(rslt)
191  class(fmsdiagoutfield_type), intent(in) :: this !< The fmsDiagOutfield_type object
192  rslt = this%phys_window
193  end function get_phys_window
194 
195  !> @brief Gets need_compute
196  !! @return copy of integer member need_compute
197  pure logical function get_need_compute (this) result(rslt)
198  class(fmsdiagoutfield_type), intent(in) :: this !< The fmsDiagOutfield_type object
199  rslt = this%need_compute
200  end function get_need_compute
201 
202  !> @brief Gets reduced_k_range
203  !! @return copy of integer member reduced_k_range
204  pure logical function get_reduced_k_range (this) result(rslt)
205  class(fmsdiagoutfield_type), intent(in) :: this !< The fmsDiagOutfield_type object
206  rslt = this%reduced_k_range
207  end function get_reduced_k_range
208 
209  !> @brief Gets missvalue_present
210  !! @return copy of integer member missvalue_present
211  pure logical function get_missvalue_present (this) result(rslt)
212  class(fmsdiagoutfield_type), intent(in) :: this !< The fmsDiagOutfield_type object
213  rslt = this%missvalue_present
214  end function get_missvalue_present
215 
216  !> @brief Gets mask_variant
217  !! @return copy of integer member mask_variant
218  pure logical function get_mask_variant (this) result(rslt)
219  class(fmsdiagoutfield_type), intent(in) :: this !< The fmsDiagOutfield_type object
220  rslt = this%mask_variant
221  end function get_mask_variant
222 
223  !> @brief Gets mask_present
224  !! @return copy of integer member mask_present
225  pure logical function get_mask_present (this) result(rslt)
226  class(fmsdiagoutfield_type), intent(in) :: this !< The fmsDiagOutfield_type object
227  rslt = this%mask_present
228  end function get_mask_present
229 
230  !> @brief Gets the time_reduction object
231  !! @return copy of the memeber object time_reduction
232  function get_time_reduction (this) result(rslt)
233  class(fmsdiagoutfield_type), intent(in) :: this !< diag object
234  TYPE(fmsdiagtimereduction_type), allocatable :: rslt
235  allocate( rslt )
236  call rslt%copy(this%time_reduction)
237  end function get_time_reduction
238 
239  !> @brief Gets f1
240  !! @return copy of integer member f1
241  pure integer function get_f1 (this) result(rslt)
242  class(fmsdiagoutfieldindex_type), intent(in) :: this !< The fmsDiagOutfieldIndex_type object
243  rslt = this%f1
244  end function get_f1
245 
246  !> @brief Gets f2
247  !! @return copy of integer member f2
248  pure integer function get_f2 (this) result(rslt)
249  class(fmsdiagoutfieldindex_type), intent(in) :: this !< The fmsDiagOutfieldIndex_type object
250  rslt = this%f2
251  end function get_f2
252 
253  !> @brief Gets f3
254  !! @return copy of integer member f3
255  pure integer function get_f3 (this) result(rslt)
256  class(fmsdiagoutfieldindex_type), intent(in) :: this !< The fmsDiagOutfieldIndex_type object
257  rslt = this%f3
258  end function get_f3
259 
260  !> @brief Gets f4
261  !! @return copy of integer member f4
262  pure integer function get_f4 (this) result(rslt)
263  class(fmsdiagoutfieldindex_type), intent(in) :: this !< The fmsDiagOutfieldIndex_type object
264  rslt = this%f4
265  end function get_f4
266 
267  !> @brief Gets is
268  !! @return copy of integer member is
269  pure integer function get_is (this) result(rslt)
270  class(fmsdiagoutfieldindex_type), intent(in) :: this !< The fmsDiagOutfieldIndex_type object
271  rslt = this%is
272  end function get_is
273 
274  !> @brief Gets js
275  !! @return copy of integer member js
276  pure integer function get_js (this) result(rslt)
277  class(fmsdiagoutfieldindex_type), intent(in) :: this !< The fmsDiagOutfieldIndex_type object
278  rslt = this%js
279  end function get_js
280 
281  !> @brief Gets ks
282  !! @return copy of integer member ks
283  pure integer function get_ks (this) result(rslt)
284  class(fmsdiagoutfieldindex_type), intent(in) :: this !< The fmsDiagOutfieldIndex_type object
285  rslt = this%ks
286  end function get_ks
287 
288  !> @brief Gets ie
289  !! @return copy of integer member ie
290  pure integer function get_ie (this) result(rslt)
291  class(fmsdiagoutfieldindex_type), intent(in) :: this !< The fmsDiagOutfieldIndex_type object
292  rslt = this%ie
293  end function get_ie
294 
295  !> @brief Gets je
296  !! @return copy of integer member je
297  pure integer function get_je (this) result(rslt)
298  class(fmsdiagoutfieldindex_type), intent(in) :: this !< The fmsDiagOutfieldIndex_type object
299  rslt = this%je
300  end function get_je
301 
302  !> @brief Gets ke
303  !! @return copy of integer member ke
304  pure integer function get_ke (this) result(rslt)
305  class(fmsdiagoutfieldindex_type), intent(in) :: this !< The fmsDiagOutfieldIndex_type object
306  rslt = this%ke
307  end function get_ke
308 
309  !> @brief Gets hi
310  !! @return copy of integer member hi
311  pure integer function get_hi (this) result(rslt)
312  class(fmsdiagoutfieldindex_type), intent(in) :: this !< The fmsDiagOutfieldIndex_type object
313  rslt = this%hi
314  end function get_hi
315 
316  !> @brief Gets hj
317  !! @return copy of integer member hj
318  pure integer function get_hj (this) result(rslt)
319  class(fmsdiagoutfieldindex_type), intent(in) :: this !< The fmsDiagOutfieldIndex_type object
320  rslt = this%hj
321  end function get_hj
322 
323 
324  !> #brief initialize all the members of the class.
325  SUBROUTINE initialize_outfield_index_type(this, is, js , ks, ie, je, ke, hi, hj, f1, f2, f3, f4)
326  CLASS(fmsdiagoutfieldindex_type), INTENT(inout) :: this
327  INTEGER, INTENT(in) :: is, js, ks !< Variable used to update class member of same names.
328  INTEGER, INTENT(in) :: ie, je, ke !< Variable used to update class member of same names.
329  INTEGER, INTENT(in) :: hi, hj !< Variable used to update class member of same names.
330  INTEGER, INTENT(in) :: f1, f2, f3, f4 !< Variable used to update class member of same names.
331 
332  this%is = is
333  this%js = js
334  this%ks = ks
335  this%ie = ie
336  this%je = je
337  this%ke = ke
338 
339  this%hi = hi
340  this%hj = hj
341 
342  this%f1 = f1
343  this%f2 = f2
344  this%f3 = f3
345  this%f4 = f4
346  END SUBROUTINE initialize_outfield_index_type
347 
348 
349  !> @brief Update the fmsDiagOutfield_type instance with those fields used in the legacy diag manager.
350  !! Note that this is initializing from the legacy structures.
351  !! Note that output_frequency came from file_type;
352  SUBROUTINE initialize_outfield_imp(this, input_field, output_field, mask_present, freq)
353  CLASS(fmsdiagoutfield_type), INTENT(inout) :: this !< An instance of the fmsDiagOutfield_type
354  TYPE(input_field_type), INTENT(in) :: input_field !< An instance of the input_field_type
355  TYPE(output_field_type), INTENT(in) :: output_field !< An instance of the output_field_type
356  LOGICAL, INTENT(in) :: mask_present !< Was the mask present in the call to send_data?
357  INTEGER, INTENT(in) :: freq !< The output frequency.
358  INTEGER :: time_redux !< The time reduction type integer.
359 
360  this%module_name = input_field%module_name
361  this%field_name = input_field%field_name
362  this%output_name = output_field%output_name
363 
364  this%pow_value = output_field%pow_value
365  this%phys_window = output_field%phys_window
366  this%need_compute = output_field%need_compute
367  this%reduced_k_range = output_field%reduced_k_range
368  this%mask_variant = input_field%mask_variant
369  !!Note: in legacy diag manager, presence of missing value vs presence of mask
370  !! is determined in different ways (diag table vs send function call)
371  this%missvalue_present = input_field%missing_value_present
372  this%mask_present = mask_present
373 
374  time_redux = get_output_field_time_reduction(output_field)
375  call this%time_reduction%initialize( time_redux , freq)
376 
377  !!TODO: the time_min and time_max buffer update code is almost the exact same src code, except
378  !! for the compariosn function. Simplify code and set comparison function:
379  !!TODO: If possible add to the power function. See issue with pointers and elemental functions
380 
381  END SUBROUTINE initialize_outfield_imp
382 
383  !> @brief Initialized an fmsDiagOutfield_type as needed for unit tests.
384  subroutine initialize_for_ut(this, module_name, field_name, output_name, &
385  & power_val, phys_window, need_compute, mask_variant, reduced_k_range, num_elems, &
386  & time_reduction_type,output_freq)
387  CLASS(fmsdiagoutfield_type), intent(inout) :: this
388  CHARACTER(len=*), INTENT(in) :: module_name !< Var with same name in fmsDiagOutfield_type
389  CHARACTER(len=*), INTENT(in) :: field_name !< Var with same name in fmsDiagOutfield_type
390  CHARACTER(len=*), INTENT(in) :: output_name !< Var with same name in fmsDiagOutfield_type
391  INTEGER, INTENT(in) :: power_val !< Var with same name in fmsDiagOutfield_type
392  LOGICAL, INTENT(in) :: phys_window !< Var with same name in fmsDiagOutfield_type
393  LOGICAL, INTENT(in) :: need_compute !< Var with same name in fmsDiagOutfield_type
394  LOGICAL, INTENT(in) :: mask_variant !< Var with same name in fmsDiagOutfield_type
395  LOGICAL, INTENT(in) :: reduced_k_range !< Var with same name in fmsDiagOutfield_type
396  INTEGER, INTENT(in) :: num_elems !< Var with same name in fmsDiagOutfield_type
397  INTEGER, INTENT(in) :: time_reduction_type !< Var with same name in fmsDiagOutfield_type
398  INTEGER, INTENT(in) :: output_freq !< The output_freq need in initaliztion of time_reduction_type
399 
400  this%module_name = module_name
401  this%field_name = field_name
402  this%output_name = output_name
403  this%pow_value = power_val
404  this%phys_window = phys_window
405  this%need_compute = need_compute
406  this%reduced_k_range = reduced_k_range
407  this%mask_variant = mask_variant
408  call this%time_reduction%initialize(time_reduction_type, output_freq)
409  end subroutine initialize_for_ut
410 
411  !> @brief Reset the time reduction member field. Intended for use in unit tests only.
412  SUBROUTINE reset_time_reduction_ut(this, source )
413  CLASS(fmsdiagoutfield_type), INTENT(inout) :: this !< An instance of the fmsDiagOutfield_type
414  TYPE(fmsdiagtimereduction_type) :: source !< The fmsDiagTimeReduction_type to copy from
415  call this%time_reduction%copy(source)
416  END SUBROUTINE reset_time_reduction_ut
417 
418 
419 
420  !> \brief Get the time reduction from a legacy output field.
421  !\note Note we do not place this in the time_reduction class to avoid circular dependencies.
422  function get_output_field_time_reduction(ofield) result (rslt)
423  TYPE(output_field_type), INTENT(in) :: ofield !< An instance of the output_field_type
424  INTEGER :: rslt !< The result integer which is the time reduction.
425  if(ofield%time_max) then
426  rslt = time_max
427  elseif(ofield%time_min)then
428  rslt = time_min
429  else if (ofield%time_sum) then
430  rslt = time_sum
431  else if (ofield%time_rms) then
432  rslt = time_rms
433  else if (ofield%time_average) then
434  rslt = time_average
435  else
436  rslt = time_none
437  !if(.NOT. ofield%static) then
438  !!TODO: Set error to FATAL. When legacy diag_manager is removed?
439  ! CALL error_mesg('fms_diag_outfield:get_output_field_time_reduction', &
440  ! & 'result is time_none but out_field%static is not true', WARNING)
441  !end if
442  endif
444 
445 END MODULE fms_diag_outfield_mod
446 !> @}
447 ! close documentation grouping
448 
449 
integer, parameter glo_reg_val_alt
Alternate value used in the region specification of the diag_table to indicate to use the full axis i...
Definition: diag_data.F90:109
logical region_out_use_alt_value
Will determine which value to use when checking a regional output if the region is the full axis or a...
Definition: diag_data.F90:377
integer, parameter glo_reg_val
Value used in the region specification of the diag_table to indicate to use the full axis instead of ...
Definition: diag_data.F90:107
Define the region for field output.
Definition: diag_data.F90:171
Type to hold the input field description.
Definition: diag_data.F90:222
Type to hold the output field description.
Definition: diag_data.F90:249
pure integer function get_ie(this)
Gets ie.
pure character(len=:) function, allocatable get_field_name(this)
Gets field_name.
integer function get_output_field_time_reduction(ofield)
Get the time reduction from a legacy output field.
pure integer function get_f1(this)
Gets f1.
subroutine initialize_outfield_imp(this, input_field, output_field, mask_present, freq)
Update the fmsDiagOutfield_type instance with those fields used in the legacy diag manager....
pure logical function get_mask_variant(this)
Gets mask_variant.
pure integer function get_js(this)
Gets js.
subroutine initialize_outfield_index_type(this, is, js, ks, ie, je, ke, hi, hj, f1, f2, f3, f4)
#brief initialize all the members of the class.
pure integer function get_ks(this)
Gets ks.
pure integer function get_f2(this)
Gets f2.
pure integer function get_hj(this)
Gets hj.
pure integer function get_je(this)
Gets je.
pure integer function get_f3(this)
Gets f3.
pure logical function get_missvalue_present(this)
Gets missvalue_present.
pure integer function get_hi(this)
Gets hi.
pure integer function get_pow_value(this)
Gets pow_value.
type(fmsdiagtimereduction_type) function, allocatable get_time_reduction(this)
Gets the time_reduction object.
pure logical function get_mask_present(this)
Gets mask_present.
subroutine reset_time_reduction_ut(this, source)
Reset the time reduction member field. Intended for use in unit tests only.
pure integer function get_ke(this)
Gets ke.
pure logical function get_need_compute(this)
Gets need_compute.
pure logical function get_reduced_k_range(this)
Gets reduced_k_range.
pure character(len=:) function, allocatable get_output_file(this)
Gets output_file.
pure integer function get_f4(this)
Gets f4.
pure character(len=:) function, allocatable get_module_name(this)
Gets module_name.
pure logical function get_phys_window(this)
Gets phys_window.
subroutine initialize_for_ut(this, module_name, field_name, output_name, power_val, phys_window, need_compute, mask_variant, reduced_k_range, num_elems, time_reduction_type, output_freq)
Initialized an fmsDiagOutfield_type as needed for unit tests.
pure character(len=:) function, allocatable get_output_name(this)
Gets output_name.
pure integer function get_is(this)
Gets is.
Class fmsDiagOutfield_type (along with class ms_diag_outfield_index_type ) contain information used i...
Class fms_diag_outfield_index_type which (along with class fmsDiagOutfield_type) encapsulate related ...
integer, parameter time_min
The reduction method is min.
integer, parameter time_power
The reduction method is power.
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.
integer, parameter time_max
The reduction method is max.
logical function, public fms_error_handler(routine, message, err_msg)
Facilitates the control of fatal error conditions.
Definition: fms.F90:468
subroutine, public error_mesg(routine, message, level)
Print notes, warnings and error messages; terminates program for warning and error messages....
Definition: fms.F90:441