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