FMS  2025.04
Flexible Modeling System
fms_diag_output_buffer.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 !> @author Ryan Mulhall
19 !> @email ryan.mulhall@noaa.gov
20 !! @brief Contains buffer types and routines for the diag manager
21 !!
22 !! @description Holds buffered data for fmsDiagVars_type objects
23 !! buffer0-5d types extend fmsDiagBuffer_class, and upon allocation
24 !! are added to the module's buffer_lists depending on it's dimension
26 #ifdef use_yaml
27 use platform_mod
28 use iso_c_binding
29 use time_manager_mod, only: time_type, operator(==), operator(>=), get_ticks_per_second, get_time, operator(>)
30 use constants_mod, only: seconds_per_day
31 use mpp_mod, only: mpp_error, fatal, note, mpp_pe, mpp_root_pe
32 use diag_data_mod, only: diag_null, diag_not_registered, i4, i8, r4, r8, get_base_time, min_value, max_value, empty, &
34 use fms2_io_mod, only: fmsnetcdffile_t, write_data, fmsnetcdfdomainfile_t, fmsnetcdfunstructureddomainfile_t
35 use fms_diag_yaml_mod, only: diag_yaml
36 use fms_diag_bbox_mod, only: fmsdiagibounds_type
37 use fms_diag_reduction_methods_mod, only: do_time_none, do_time_min, do_time_max, do_time_sum_update, time_update_done
38 use fms_diag_time_utils_mod, only: diag_time_inc
39 
40 implicit none
41 
42 private
43 
44 !> holds an allocated buffer0-5d object
46  integer :: buffer_id !< index in buffer list
47  integer(i4_kind) :: buffer_type !< set to allocated data type & kind value, one of i4,i8,r4,r8
48  class(*), allocatable :: buffer(:,:,:,:,:) !< 5D numeric data array
49  integer :: ndim !< Number of dimensions for each variable
50  integer, allocatable :: buffer_dims(:) !< holds the size of each dimension in the buffer
51  real(r8_kind), allocatable :: weight_sum(:,:,:,:) !< Weight sum as an array
52  !! (this will be have a size of 1,1,1,1 when not using variable
53  !! masks!)
54  integer, allocatable :: num_elements(:) !< used in time-averaging
55  integer, allocatable :: axis_ids(:) !< Axis ids for the buffer
56  integer :: field_id !< The id of the field the buffer belongs to
57  integer :: yaml_id !< The id of the yaml id the buffer belongs to
58  logical :: done_with_math !< .True. if done doing the math
59  integer :: diurnal_sample_size = -1 !< dirunal sample size as read in from the reduction method
60  !! ie. diurnal24 = sample size of 24
61  integer :: diurnal_section= -1 !< the diurnal section (ie 5th index) calculated from the current model
62  !! time and sample size if using a diurnal reduction
63  logical, allocatable :: send_data_called !< .True. if send_data has been called
64  integer :: unlmited_dimension !< Unlimited dimension index of the last write for this output buffer
65  type(time_type) :: time !< The last time the data was received
66  type(time_type) :: next_output !< The next time to output the data
67 
68  contains
69  procedure :: get_buffer_name
70  procedure :: add_axis_ids
71  procedure :: get_axis_ids
72  procedure :: set_field_id
73  procedure :: get_field_id
74  procedure :: set_yaml_id
75  procedure :: get_yaml_id
76  procedure :: init_buffer_time
77  procedure :: set_next_output
78  procedure :: update_buffer_time
79  procedure :: get_buffer_time
80  procedure :: is_there_data_to_write
81  procedure :: is_time_to_finish_reduction
82  procedure :: set_send_data_called
83  procedure :: is_done_with_math
84  procedure :: set_done_with_math
85  procedure :: write_buffer
86  procedure :: init_buffer_unlim_dim
87  procedure :: increase_unlim_dim
88  procedure :: get_unlim_dim
89  !! These are needed because otherwise the write_data calls will go into the wrong interface
90  procedure :: write_buffer_wrapper_netcdf
91  procedure :: write_buffer_wrapper_domain
92  procedure :: write_buffer_wrapper_u
93  procedure :: allocate_buffer
94  procedure :: initialize_buffer
95  procedure :: get_buffer
96  procedure :: flush_buffer
97  procedure :: do_time_none_wrapper
98  procedure :: do_time_min_wrapper
99  procedure :: do_time_max_wrapper
100  procedure :: do_time_sum_wrapper
101  procedure :: diag_reduction_done_wrapper
102  procedure :: get_buffer_dims
103  procedure :: get_diurnal_sample_size
104  procedure :: set_diurnal_sample_size
105  procedure :: set_diurnal_section_index
106  procedure :: get_remapped_diurnal_data
108 
109 ! public types
110 public :: fmsdiagoutputbuffer_type
111 
112 ! public routines
113 public :: fms_diag_output_buffer_init
114 
115 contains
116 
117 !!--------module routines
118 
119 !> Initializes a list of diag buffers
120 !> @returns true if allocation is successfull
121 logical function fms_diag_output_buffer_init(buffobjs, buff_list_size)
122  type(fmsdiagoutputbuffer_type), allocatable, intent(out) :: buffobjs(:) !< an array of buffer container types
123  !! to allocate
124  integer, intent(in) :: buff_list_size !< size of buffer array to allocate
125 
126  if (allocated(buffobjs)) call mpp_error(fatal,'fms_diag_buffer_init: passed in buffobjs array is already allocated')
127  allocate(buffobjs(buff_list_size))
128  fms_diag_output_buffer_init = allocated(buffobjs)
129 end function fms_diag_output_buffer_init
130 
131 !!--------generic routines for any fmsDiagBuffer_class objects
132 
133 !> Setter for buffer_id for any buffer objects
134 subroutine set_buffer_id(this, id)
135  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< buffer object to set id for
136  integer, intent(in) :: id !< positive integer id to set
137 
138  this%buffer_id = id
139 end subroutine set_buffer_id
140 
141 !> Deallocates data fields from a buffer object.
142 subroutine flush_buffer(this)
143  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< any buffer object
144 
145  this%buffer_id = diag_null
146  this%buffer_type = diag_null
147  this%ndim = diag_null
148  this%field_id = diag_null
149  this%yaml_id = diag_null
150  if (allocated(this%buffer)) deallocate(this%buffer)
151  if (allocated(this%buffer_dims)) deallocate(this%buffer_dims)
152  if (allocated(this%num_elements)) deallocate(this%num_elements)
153  if (allocated(this%axis_ids)) deallocate(this%axis_ids)
154  if (allocated(this%weight_sum)) deallocate(this%weight_sum)
155 end subroutine flush_buffer
156 
157 !> Allocates a 5D buffer to given buff_type.
158 subroutine allocate_buffer(this, buff_type, ndim, buff_sizes, mask_variant, field_name, diurnal_samples)
159  class(fmsdiagoutputbuffer_type), intent(inout), target :: this !< 5D buffer object
160  class(*), intent(in) :: buff_type !< allocates to the type of buff_type
161  integer, intent(in) :: ndim !< Number of dimension
162  integer, intent(in) :: buff_sizes(4) !< dimension buff_sizes
163  logical, intent(in) :: mask_variant !< Mask changes over time
164  character(len=*), intent(in) :: field_name !< field name for error output
165  integer, intent(in) :: diurnal_samples !< number of diurnal samples
166 
167  integer :: n_samples !< number of diurnal samples, defaults to 1
168 
169  n_samples = max(1, diurnal_samples)
170  call this%set_diurnal_sample_size(n_samples)
171 
172  this%ndim =ndim
173  if(allocated(this%buffer)) call mpp_error(fatal, "allocate_buffer: buffer already allocated for field:" // &
174  field_name)
175  select type (buff_type)
176  type is (integer(kind=i4_kind))
177  allocate(integer(kind=i4_kind) :: this%buffer(buff_sizes(1),buff_sizes(2),buff_sizes(3),buff_sizes(4), &
178  & n_samples))
179  this%buffer_type = i4
180  type is (integer(kind=i8_kind))
181  allocate(integer(kind=i8_kind) :: this%buffer(buff_sizes(1),buff_sizes(2),buff_sizes(3),buff_sizes(4), &
182  & n_samples))
183  this%buffer_type = i8
184  type is (real(kind=r4_kind))
185  allocate(real(kind=r4_kind) :: this%buffer(buff_sizes(1),buff_sizes(2),buff_sizes(3),buff_sizes(4), &
186  & n_samples))
187  this%buffer_type = r4
188  type is (real(kind=r8_kind))
189  allocate(real(kind=r8_kind) :: this%buffer(buff_sizes(1),buff_sizes(2),buff_sizes(3),buff_sizes(4), &
190  & n_samples))
191  this%buffer_type = r8
192  class default
193  call mpp_error("allocate_buffer", &
194  "The buff_type value passed to allocate a buffer is not a r8, r4, i8, or i4" // &
195  "for field:" // field_name, fatal)
196  end select
197  if (mask_variant) then
198  allocate(this%weight_sum(buff_sizes(1),buff_sizes(2),buff_sizes(3),buff_sizes(4)))
199  else
200  allocate(this%weight_sum(1,1,1,1))
201  endif
202  this%weight_sum = 0.0_r8_kind
203 
204  allocate(this%num_elements(n_samples))
205  this%num_elements = 0
206  this%done_with_math = .false.
207  this%send_data_called = .false.
208  allocate(this%buffer_dims(5))
209  this%buffer_dims(1) = buff_sizes(1)
210  this%buffer_dims(2) = buff_sizes(2)
211  this%buffer_dims(3) = buff_sizes(3)
212  this%buffer_dims(4) = buff_sizes(4)
213  this%buffer_dims(5) = n_samples
214 end subroutine allocate_buffer
215 
216 !> Get routine for 5D buffers.
217 !! Sets the buff_out argument to the integer or real array currently stored in the buffer.
218 subroutine get_buffer (this, buff_out, field_name)
219  class(fmsdiagoutputbuffer_type), intent(in) :: this !< 5d allocated buffer object
220  class(*), allocatable, intent(out) :: buff_out(:,:,:,:,:) !< output of copied buffer data
221  !! must be the same size as the allocated buffer
222  character(len=*), intent(in) :: field_name !< field name for error output
223 
224  integer(i4_kind) :: buff_size(5)!< sizes for allocated buffer
225 
226  if(.not. allocated(this%buffer)) call mpp_error(fatal, 'get_buffer: buffer not yet allocated for field:' &
227  & // field_name)
228  buff_size(1) = size(this%buffer,1)
229  buff_size(2) = size(this%buffer,2)
230  buff_size(3) = size(this%buffer,3)
231  buff_size(4) = size(this%buffer,4)
232  buff_size(5) = size(this%buffer,5)
233 
234  select type (buff=>this%buffer)
235  type is (real(r4_kind))
236  allocate(real(r4_kind) :: buff_out(buff_size(1), buff_size(2), buff_size(3), buff_size(4), buff_size(5)))
237  buff_out = buff
238  type is (real(r8_kind))
239  allocate(real(r8_kind) :: buff_out(buff_size(1), buff_size(2), buff_size(3), buff_size(4), buff_size(5)))
240  buff_out = buff
241  type is (integer(i4_kind))
242  allocate(integer(i4_kind) :: buff_out(buff_size(1), buff_size(2), buff_size(3), buff_size(4), buff_size(5)))
243  buff_out = buff
244  type is (integer(i8_kind))
245  allocate(integer(i8_kind) :: buff_out(buff_size(1), buff_size(2), buff_size(3), buff_size(4), buff_size(5)))
246  buff_out = buff
247  class default
248  call mpp_error(fatal, "get_buffer: buffer allocated to invalid type(must be integer or real, kind size 4 or 8)."&
249  //"field name: "// field_name)
250  end select
251 end subroutine
252 
253 !> @brief Initializes a buffer based on the reduction method
254 subroutine initialize_buffer (this, reduction_method, field_name)
255  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< allocated 5D buffer object
256  integer, intent(in) :: reduction_method !< The reduction method for the field
257  character(len=*), intent(in) :: field_name !< field name for error output
258 
259  if(.not. allocated(this%buffer)) call mpp_error(fatal, 'initialize_buffer: field:'// field_name // &
260  'buffer not yet allocated, allocate_buffer() must be called on this object first.')
261 
262  select type(buff => this%buffer)
263  type is(real(r8_kind))
264  select case (reduction_method)
265  case (time_min)
266  buff = real(min_value, kind=r8_kind)
267  case (time_max)
268  buff = real(max_value, kind=r8_kind)
269  case default
270  buff = real(empty, kind=r8_kind)
271  end select
272  type is(real(r4_kind))
273  select case (reduction_method)
274  case (time_min)
275  buff = real(min_value, kind=r4_kind)
276  case (time_max)
277  buff = real(max_value, kind=r4_kind)
278  case default
279  buff = real(empty, kind=r4_kind)
280  end select
281  type is(integer(i8_kind))
282  select case (reduction_method)
283  case (time_min)
284  buff = int(min_value, kind=i8_kind)
285  case (time_max)
286  buff = int(max_value, kind=i8_kind)
287  case default
288  buff = int(empty, kind=i8_kind)
289  end select
290  type is(integer(i4_kind))
291  select case (reduction_method)
292  case (time_min)
293  buff = int(min_value, kind=i4_kind)
294  case (time_max)
295  buff = int(max_value, kind=i4_kind)
296  case default
297  buff = int(empty, kind=i4_kind)
298  end select
299  class default
300  call mpp_error(fatal, 'initialize buffer_5d: buffer allocated to invalid data type, this shouldnt happen')
301  end select
302 
303 end subroutine initialize_buffer
304 
305 !> @brief Get the name of the field for the output buffer
306 !! @return Name of the field for the output buffer
307 function get_buffer_name(this) &
308  result(rslt)
309  class(fmsdiagoutputbuffer_type), intent(in) :: this !< Buffer object
310 
311  character(len=:), allocatable :: rslt
312 
313  rslt = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
314 end function get_buffer_name
315 
316 !> @brief Adds the axis ids to the buffer object
317 subroutine add_axis_ids(this, axis_ids)
318  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< Buffer object
319  integer, intent(in) :: axis_ids(:) !< Axis ids to add
320 
321  this%axis_ids = axis_ids
322 end subroutine
323 
324 !> @brief Get the axis_ids for the buffer
325 !! @return Axis_ids, if the buffer doesn't have axis ids it returns diag_null
326 subroutine get_axis_ids(this, res)
327  class(fmsdiagoutputbuffer_type), target, intent(inout) :: this !< Buffer object
328  integer, pointer, intent(out) :: res(:)
329 
330  if (allocated(this%axis_ids)) then
331  res => this%axis_ids
332  else
333  allocate(res(1))
334  res = diag_null
335  endif
336 end subroutine
337 
338 !> @brief Get the field id of the buffer
339 !! @return the field id of the buffer
340 function get_field_id(this) &
341  result(res)
342  class(fmsdiagoutputbuffer_type), intent(in) :: this !< Buffer object
343  integer :: res
344 
345  res = this%field_id
346 end function get_field_id
347 
348 !> @brief set the field id of the buffer
349 subroutine set_field_id(this, field_id)
350  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< Buffer object
351  integer, intent(in) :: field_id !< field id of the buffer
352 
353  this%field_id = field_id
354 end subroutine set_field_id
355 
356 !> @brief set the field id of the buffer
357 subroutine set_yaml_id(this, yaml_id)
358  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< Buffer object
359  integer, intent(in) :: yaml_id !< yaml id of the buffer
360 
361  this%yaml_id = yaml_id
362 end subroutine set_yaml_id
363 
364 !> @brief inits the buffer time for the buffer
365 subroutine init_buffer_time(this, time)
366  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< Buffer object
367  type(time_type), optional, intent(in) :: time !< time to add to the buffer
368 
369  if (present(time)) then
370  this%time = time
371  this%next_output = time
372  else
373  this%time = get_base_time()
374  this%next_output = this%time
375  endif
376 end subroutine init_buffer_time
377 
378 !> @brief Sets the next output
379 subroutine set_next_output(this, next_output, next_next_output, is_static)
380  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< Buffer object
381  type(time_type), intent(in) :: next_output !< The current next_output in the file obj
382  type(time_type), intent(in) :: next_next_output !< The current next_next_output in the file obj
383  logical, optional, intent(in) :: is_static !< .True. if the field is static
384 
385  if (present(is_static)) then
386  !< If the field is static set the next_output to be equal to time
387  !! this should only be used in the init, so next_output will be equal to the init time
388  if (is_static) then
389  this%next_output = this%time
390  return
391  endif
392  endif
393 
394  !< If the file's next_output is greater than the buffer's next output set
395  !! the buffer's next output to the file's next_ouput, otherwise use the file's
396  !! next_next_output
397  !! This is needed for when file have fields that get data send data sent at different frequencies
398  if (next_output > this%next_output) then
399  this%next_output = next_output
400  else
401  this%next_output = next_next_output
402  endif
403 end subroutine set_next_output
404 
405 !> @brief Update the buffer time if it is a new time
406 subroutine update_buffer_time(this, time)
407  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< Buffer object
408  type(time_type), intent(in) :: time !< Current model time
409 
410  if (time > this%time) then
411  this%time = time
412  endif
413 end subroutine update_buffer_time
414 
415 !> @brief Get the buffer_time from a output buffer object
416 !! @return The buffer time
417 function get_buffer_time(this) &
418  result(rslt)
419  class(fmsdiagoutputbuffer_type), intent(in) :: this !< Buffer object
420  type(time_type) :: rslt
421 
422  rslt = this%time
423 end function get_buffer_time
424 
425 !> @brief Determine if finished with math
426 !! @return this%done_with_math
427 function is_done_with_math(this) &
428  result(res)
429  class(fmsdiagoutputbuffer_type), intent(in) :: this !< Buffer object
430  logical :: res
431 
432  res = this%done_with_math
433 end function is_done_with_math
434 
435 !> @brief Set done_with_math to .true.
436 subroutine set_done_with_math(this)
437  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< Buffer object
438  integer :: res
439 
440  this%done_with_math = .true.
441 end subroutine set_done_with_math
442 
443 !> @brief Get the yaml id of the buffer
444 !! @return the yaml id of the buffer
445 function get_yaml_id(this) &
446  result(res)
447 
448  class(fmsdiagoutputbuffer_type), intent(in) :: this !< Buffer object
449  integer :: res
450 
451  res = this%yaml_id
452 end function get_yaml_id
453 
454 !> @brief Get the unlim dimension index of the buffer object
455 !! @return The unlim dimension index of the buffer object
456 function get_unlim_dim(this) &
457  result(res)
458  class(fmsdiagoutputbuffer_type), intent(in) :: this !< buffer object to write
459  integer :: res
460 
461  res = this%unlmited_dimension
462 end function get_unlim_dim
463 
464 !> @brief Increase the unlim dimension index of the buffer object
465 subroutine increase_unlim_dim(this)
466  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< buffer object to write
467 
468  this%unlmited_dimension = this%unlmited_dimension + 1
469 end subroutine increase_unlim_dim
470 
471 !> @brief Init the unlim dimension index of the buffer object to 0
472 subroutine init_buffer_unlim_dim(this)
473  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< buffer object to write
474 
475  this%unlmited_dimension = 0
476 end subroutine init_buffer_unlim_dim
477 
478 !> @brief Write the buffer to the file
479 subroutine write_buffer(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
480  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< buffer object to write
481  class(fmsnetcdffile_t), intent(in) :: fms2io_fileobj !< fileobj to write to
482  integer, optional, intent(in) :: unlim_dim_level !< unlimited dimension
483  logical, optional, intent(in) :: is_diurnal !< should be set if using diurnal
484  !! reductions so buffer data can be remapped
485 
486  select type(fms2io_fileobj)
487  type is (fmsnetcdffile_t)
488  call this%write_buffer_wrapper_netcdf(fms2io_fileobj, unlim_dim_level=unlim_dim_level, is_diurnal=is_diurnal)
489  type is (fmsnetcdfdomainfile_t)
490  call this%write_buffer_wrapper_domain(fms2io_fileobj, unlim_dim_level=unlim_dim_level, is_diurnal=is_diurnal)
491  type is (fmsnetcdfunstructureddomainfile_t)
492  call this%write_buffer_wrapper_u(fms2io_fileobj, unlim_dim_level=unlim_dim_level, is_diurnal=is_diurnal)
493  class default
494  call mpp_error(fatal, "The file "//trim(fms2io_fileobj%path)//" is not one of the accepted types"//&
495  " only FmsNetcdfFile_t, FmsNetcdfDomainFile_t, and FmsNetcdfUnstructuredDomainFile_t are accepted.")
496  end select
497 
498  call this%initialize_buffer(diag_yaml%diag_fields(this%yaml_id)%get_var_reduction(), &
499  diag_yaml%diag_fields(this%yaml_id)%get_var_outname())
500  !TODO Set the counters back to 0
501 end subroutine write_buffer
502 
503 !> @brief Write the buffer to the FmsNetcdfFile_t fms2io_fileobj
504 subroutine write_buffer_wrapper_netcdf(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
505  class(fmsdiagoutputbuffer_type), intent(in) :: this !< buffer object to write
506  type(fmsnetcdffile_t), intent(in) :: fms2io_fileobj !< fileobj to write to
507  integer, optional, intent(in) :: unlim_dim_level !< unlimited dimension
508  logical, optional, intent(in) :: is_diurnal !< should be set if using diurnal
509  !! reductions so buffer data can be remapped
510  character(len=:), allocatable :: varname !< name of the variable
511  logical :: using_diurnal !< local copy of is_diurnal if present
512  class(*), allocatable :: buff_ptr(:,:,:,:,:) !< pointer for buffer to write
513 
514  using_diurnal = .false.
515  if( present(is_diurnal) ) using_diurnal = is_diurnal
516  if( using_diurnal ) then
517  call this%get_remapped_diurnal_data(buff_ptr)
518  else
519  buff_ptr = this%buffer
520  endif
521 
522  varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
523  select case(this%ndim)
524  case (0)
525  call write_data(fms2io_fileobj, varname, buff_ptr(1,1,1,1,1), unlim_dim_level=unlim_dim_level)
526  case (1)
527  call write_data(fms2io_fileobj, varname, buff_ptr(:,1,1,1,1), unlim_dim_level=unlim_dim_level)
528  case (2)
529  call write_data(fms2io_fileobj, varname, buff_ptr(:,:,1,1,1), unlim_dim_level=unlim_dim_level)
530  case (3)
531  call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,1,1), unlim_dim_level=unlim_dim_level)
532  case (4)
533  call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,1), unlim_dim_level=unlim_dim_level)
534  case (5)
535  call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,:), unlim_dim_level=unlim_dim_level)
536  end select
537 end subroutine write_buffer_wrapper_netcdf
538 
539 !> @brief Write the buffer to the FmsNetcdfDomainFile_t fms2io_fileobj
540 subroutine write_buffer_wrapper_domain(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
541  class(fmsdiagoutputbuffer_type), intent(in) :: this !< buffer object to write
542  type(fmsnetcdfdomainfile_t), intent(in) :: fms2io_fileobj !< fileobj to write to
543  integer, optional, intent(in) :: unlim_dim_level !< unlimited dimension
544  logical, optional, intent(in) :: is_diurnal !< should be set if using diurnal
545  !! reductions so buffer data can be remapped
546 
547  character(len=:), allocatable :: varname !< name of the variable
548  logical :: using_diurnal !< local copy of is_diurnal if present
549  class(*), allocatable :: buff_ptr(:,:,:,:,:) !< pointer to buffer to write
550 
551  using_diurnal = .false.
552  if( present(is_diurnal) ) using_diurnal = is_diurnal
553  if( using_diurnal ) then
554  call this%get_remapped_diurnal_data(buff_ptr)
555  else
556  buff_ptr = this%buffer
557  endif
558 
559  varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
560  select case(this%ndim)
561  case (0)
562  call write_data(fms2io_fileobj, varname, buff_ptr(1,1,1,1,1), unlim_dim_level=unlim_dim_level)
563  case (1)
564  call write_data(fms2io_fileobj, varname, buff_ptr(:,1,1,1,1), unlim_dim_level=unlim_dim_level)
565  case (2)
566  call write_data(fms2io_fileobj, varname, buff_ptr(:,:,1,1,1), unlim_dim_level=unlim_dim_level)
567  case (3)
568  call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,1,1), unlim_dim_level=unlim_dim_level)
569  case (4)
570  call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,1), unlim_dim_level=unlim_dim_level)
571  case (5)
572  call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,:), unlim_dim_level=unlim_dim_level)
573  end select
574 end subroutine write_buffer_wrapper_domain
575 
576 !> @brief Write the buffer to the FmsNetcdfUnstructuredDomainFile_t fms2io_fileobj
577 subroutine write_buffer_wrapper_u(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
578  class(fmsdiagoutputbuffer_type), intent(in) :: this !< buffer object to write
579  type(fmsnetcdfunstructureddomainfile_t), intent(in) :: fms2io_fileobj !< fileobj to write to
580  integer, optional, intent(in) :: unlim_dim_level !< unlimited dimension
581  logical, optional, intent(in) :: is_diurnal !< should be set if using diurnal
582  !! reductions so buffer data can be remapped
583 
584  character(len=:), allocatable :: varname !< name of the variable
585  logical :: using_diurnal !< local copy of is_diurnal if present
586  class(*), allocatable :: buff_ptr(:,:,:,:,:) !< pointer for buffer to write
587 
588  using_diurnal = .false.
589  if( present(is_diurnal) ) using_diurnal = is_diurnal
590  if( using_diurnal ) then
591  call this%get_remapped_diurnal_data(buff_ptr)
592  else
593  buff_ptr = this%buffer
594  endif
595 
596  varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
597  select case(this%ndim)
598  case (0)
599  call write_data(fms2io_fileobj, varname, buff_ptr(1,1,1,1,1), unlim_dim_level=unlim_dim_level)
600  case (1)
601  call write_data(fms2io_fileobj, varname, buff_ptr(:,1,1,1,1), unlim_dim_level=unlim_dim_level)
602  case (2)
603  call write_data(fms2io_fileobj, varname, buff_ptr(:,:,1,1,1), unlim_dim_level=unlim_dim_level)
604  case (3)
605  call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,1,1), unlim_dim_level=unlim_dim_level)
606  case (4)
607  call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,1), unlim_dim_level=unlim_dim_level)
608  case (5)
609  call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,:), unlim_dim_level=unlim_dim_level)
610  end select
611 end subroutine write_buffer_wrapper_u
612 
613 !> @brief Does the time_none reduction method on the buffer object
614 !! @return Error message if the math was not successful
615 function do_time_none_wrapper(this, field_data, mask, is_masked, bounds_in, bounds_out, missing_value) &
616  result(err_msg)
617  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< buffer object to write
618  class(*), intent(in) :: field_data(:,:,:,:) !< Buffer data for current time
619  type(fmsdiagibounds_type), intent(in) :: bounds_in !< Indicies for the buffer passed in
620  type(fmsdiagibounds_type), intent(in) :: bounds_out !< Indicies for the output buffer
621  logical, intent(in) :: mask(:,:,:,:) !< Mask for the field
622  logical, intent(in) :: is_masked !< .True. if the field has a mask
623  real(kind=r8_kind), intent(in) :: missing_value !< Missing_value for data points that are masked
624  character(len=50) :: err_msg
625 
626  !TODO This will be expanded for integers
627  err_msg = ""
628  select type (output_buffer => this%buffer)
629  type is (real(kind=r8_kind))
630  select type (field_data)
631  type is (real(kind=r8_kind))
632  call do_time_none(output_buffer, field_data, mask, is_masked, bounds_in, bounds_out, missing_value)
633  class default
634  err_msg="do_time_none_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
635  end select
636  type is (real(kind=r4_kind))
637  select type (field_data)
638  type is (real(kind=r4_kind))
639  call do_time_none(output_buffer, field_data, mask, is_masked, bounds_in, bounds_out, &
640  real(missing_value, kind=r4_kind))
641  class default
642  err_msg="do_time_none_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
643  end select
644  end select
645 end function do_time_none_wrapper
646 
647 !> @brief Does the time_min reduction method on the buffer object
648 !! @return Error message if the math was not successful
649 function do_time_min_wrapper(this, field_data, mask, is_masked, bounds_in, bounds_out, missing_value) &
650  result(err_msg)
651  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< buffer object to write
652  class(*), intent(in) :: field_data(:,:,:,:) !< Buffer data for current time
653  type(fmsdiagibounds_type), intent(in) :: bounds_in !< Indicies for the buffer passed in
654  type(fmsdiagibounds_type), intent(in) :: bounds_out !< Indicies for the output buffer
655  logical, intent(in) :: mask(:,:,:,:) !< Mask for the field
656  logical, intent(in) :: is_masked !< .True. if the field has a mask
657  real(kind=r8_kind), intent(in) :: missing_value !< Missing_value for data points that are masked
658  character(len=50) :: err_msg
659 
660  !TODO This will be expanded for integers
661  err_msg = ""
662  select type (output_buffer => this%buffer)
663  type is (real(kind=r8_kind))
664  select type (field_data)
665  type is (real(kind=r8_kind))
666  call do_time_min(output_buffer, field_data, mask, is_masked, bounds_in, bounds_out, missing_value)
667  class default
668  err_msg="do_time_min_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
669  end select
670  type is (real(kind=r4_kind))
671  select type (field_data)
672  type is (real(kind=r4_kind))
673  call do_time_min(output_buffer, field_data, mask, is_masked, bounds_in, bounds_out, &
674  real(missing_value, kind=r4_kind))
675  class default
676  err_msg="do_time_min_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
677  end select
678  end select
679 end function do_time_min_wrapper
680 
681 !> @brief Does the time_min reduction method on the buffer object
682 !! @return Error message if the math was not successful
683 function do_time_max_wrapper(this, field_data, mask, is_masked, bounds_in, bounds_out, missing_value) &
684  result(err_msg)
685  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< buffer object to write
686  class(*), intent(in) :: field_data(:,:,:,:) !< Buffer data for current time
687  type(fmsdiagibounds_type), intent(in) :: bounds_in !< Indicies for the buffer passed in
688  type(fmsdiagibounds_type), intent(in) :: bounds_out !< Indicies for the output buffer
689  logical, intent(in) :: mask(:,:,:,:) !< Mask for the field
690  logical, intent(in) :: is_masked !< .True. if the field has a mask
691  real(kind=r8_kind), intent(in) :: missing_value !< Missing_value for data points that are masked
692  character(len=50) :: err_msg
693 
694  !TODO This will be expanded for integers
695  err_msg = ""
696  select type (output_buffer => this%buffer)
697  type is (real(kind=r8_kind))
698  select type (field_data)
699  type is (real(kind=r8_kind))
700  call do_time_max(output_buffer, field_data, mask, is_masked, bounds_in, bounds_out, missing_value)
701  class default
702  err_msg="do_time_max_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
703  end select
704  type is (real(kind=r4_kind))
705  select type (field_data)
706  type is (real(kind=r4_kind))
707  call do_time_max(output_buffer, field_data, mask, is_masked, bounds_in, bounds_out, &
708  real(missing_value, kind=r4_kind))
709  class default
710  err_msg="do_time_max_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
711  end select
712  end select
713 end function do_time_max_wrapper
714 
715 !> @brief Does the time_sum reduction method on the buffer object
716 !! @return Error message if the math was not successful
717 function do_time_sum_wrapper(this, field_data, mask, is_masked, mask_variant, bounds_in, bounds_out, missing_value, &
718  has_missing_value, pow_value, weight) &
719  result(err_msg)
720  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< buffer object to write
721  class(*), intent(in) :: field_data(:,:,:,:) !< Buffer data for current time
722  type(fmsdiagibounds_type), intent(in) :: bounds_in !< Indicies for the buffer passed in
723  type(fmsdiagibounds_type), intent(in) :: bounds_out !< Indicies for the output buffer
724  logical, intent(in) :: mask(:,:,:,:) !< Mask for the field
725  logical, intent(in) :: is_masked !< .True. if the field has a mask
726  logical, intent(in) :: mask_variant !< .True. if the mask changes over time
727  real(kind=r8_kind), intent(in) :: missing_value !< Missing_value for data points that are masked
728  logical, intent(in) :: has_missing_value !< .True. if the field was registered with
729  !! a missing value
730  integer, optional, intent(in) :: pow_value !< power value, will calculate field_data^pow
731  !! before adding to buffer should only be
732  !! present if using pow reduction method
733  real(kind=r8_kind), optional, intent(in) :: weight !< The weight to use when suming
734  character(len=150) :: err_msg
735 
736  !TODO This will be expanded for integers
737  err_msg = ""
738  select type (output_buffer => this%buffer)
739  type is (real(kind=r8_kind))
740  select type (field_data)
741  type is (real(kind=r8_kind))
742  if (.not. is_masked) then
743  if (any(field_data .eq. missing_value)) &
744  err_msg = "You cannot pass data with missing values without masking them!"
745  endif
746  call do_time_sum_update(output_buffer, this%weight_sum, field_data, mask, is_masked, mask_variant, &
747  bounds_in, bounds_out, missing_value, this%diurnal_section, &
748  pow=pow_value, weight=weight)
749  class default
750  err_msg="do_time_sum_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
751  end select
752  type is (real(kind=r4_kind))
753  select type (field_data)
754  type is (real(kind=r4_kind))
755  if (.not. is_masked) then
756  if (any(field_data .eq. missing_value)) &
757  err_msg = "You cannot pass data with missing values without masking them!"
758  endif
759  call do_time_sum_update(output_buffer, this%weight_sum, field_data, mask, is_masked, mask_variant, &
760  bounds_in, bounds_out, real(missing_value, kind=r4_kind), &
761  this%diurnal_section, pow=pow_value, weight=weight)
762  class default
763  err_msg="do_time_sum_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
764  end select
765  class default
766  err_msg="do_time_sum_wrapper::the output buffer is not a valid type, must be real(r8_kind) or real(r4_kind)"
767  end select
768 end function do_time_sum_wrapper
769 
770 !> Finishes calculations for any reductions that use an average (avg, rms, pow)
771 !! TODO add mask and any other needed args for adjustment, and pass in the adjusted mask
772 !! to time_update_done
773 function diag_reduction_done_wrapper(this, reduction_method, missing_value, has_mask, mask_variant) &
774  result(err_msg)
775  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< Updated buffer object
776  integer, intent(in) :: reduction_method !< enumerated reduction type from diag_data
777  real(kind=r8_kind), intent(in) :: missing_value !< missing_value for masked data points
778  logical, intent(in) :: has_mask !< indicates if there was a mask used during buffer updates
779  logical, intent(in) :: mask_variant !< Indicates if the mask changes over time
780  character(len=51) :: err_msg !< error message to return, blank if sucessful
781 
782  if(.not. allocated(this%buffer)) return
783 
784  err_msg = ""
785  select type(buff => this%buffer)
786  type is (real(r8_kind))
787  call time_update_done(buff, this%weight_sum, reduction_method, missing_value, has_mask, mask_variant, &
788  this%diurnal_sample_size)
789  type is (real(r4_kind))
790  call time_update_done(buff, this%weight_sum, reduction_method, real(missing_value, r4_kind), has_mask, &
791  mask_variant, this%diurnal_sample_size)
792  end select
793  this%weight_sum = 0.0_r8_kind
794 
795 end function
796 
797 !> this leaves out the diurnal index cause its only used for tmp mask allocation
798 pure function get_buffer_dims(this)
799  class(fmsdiagoutputbuffer_type), intent(in) :: this !< buffer object to get from
800  integer :: get_buffer_dims(4)
801  get_buffer_dims = this%buffer_dims(1:4)
802 end function
803 
804 !> Get diurnal sample size (amount of diurnal sections)
805 pure integer function get_diurnal_sample_size(this)
806  class(fmsdiagoutputbuffer_type), intent(in) :: this !< buffer object to get from
807  get_diurnal_sample_size = this%diurnal_sample_size
808 end function get_diurnal_sample_size
809 
810 !> Set diurnal sample size (amount of diurnal sections)
811 subroutine set_diurnal_sample_size(this, sample_size)
812  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< buffer object to set sample size for
813  integer, intent(in) :: sample_size !< sample size to used to split daily
814  !! data into given amount of sections
815  this%diurnal_sample_size = sample_size
816 end subroutine set_diurnal_sample_size
817 
818 !> Set diurnal section index based off the current time and previously set diurnal_samplesize
819 !! Calculates which diurnal section of daily data the current time is in
820 subroutine set_diurnal_section_index(this, time)
821  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< buffer object to set diurnal index for
822  type(time_type), intent(in) :: time !< current model time
823  integer :: seconds, days, ticks
824 
825  if(this%diurnal_sample_size .lt. 0) call mpp_error(fatal, "set_diurnal_section_index::"// &
826  " diurnal sample size must be set before trying to set diurnal index for send_data")
827 
828  call get_time(time,seconds,days,ticks) ! get current date
829  ! calculates which diurnal section current time is in for a given amount of diurnal sections(<24)
830  this%diurnal_section = floor( (seconds+real(ticks)/get_ticks_per_second()) &
831  & * this%diurnal_sample_size/seconds_per_day) + 1
832 end subroutine set_diurnal_section_index
833 
834 !> Remaps the output buffer array when using the diurnal reduction
835 !! moves the diurnal index to the left-most unused dimension for the io
836 subroutine get_remapped_diurnal_data(this, res)
837  class(fmsdiagoutputbuffer_type), intent(in) :: this !< output buffer object
838  class(*), intent(out), allocatable :: res(:,:,:,:,:) !< resulting remapped data
839  integer :: last_dim !< last dimension thats used
840  integer :: ie, je, ke, ze, de !< ending indices for the new array
841  integer(i4_kind) :: buff_size(5)!< sizes for allocated buffer
842 
843  ! last dim is number of dimensions - 1 for diurnal axis
844  last_dim = this%ndim - 1
845  ! get the bounds of the remapped output array based on # of dims
846  ke = 1; ze = 1; de = 1
847  select case(last_dim)
848  case (1)
849  ie = this%buffer_dims(1); je = this%buffer_dims(5)
850  case (2)
851  ie = this%buffer_dims(1); je = this%buffer_dims(2)
852  ke = this%buffer_dims(5)
853  case (3)
854  ie = this%buffer_dims(1); je = this%buffer_dims(2)
855  ke = this%buffer_dims(3); ze = this%buffer_dims(5)
856  case (4)
857  ! no need to remap if 4d
858  res = this%buffer
859  return
860  end select
861 
862  select type(buff => this%buffer)
863  type is (real(r8_kind))
864  allocate(real(r8_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
865  select type(res)
866  type is (real(r8_kind))
867  res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, shape(res))
868  end select
869  type is (real(r4_kind))
870  allocate(real(r4_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
871  select type(res)
872  type is (real(r4_kind))
873  res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, shape(res))
874  end select
875  type is (integer(i8_kind))
876  allocate(integer(i8_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
877  select type(res)
878  type is (integer(i8_kind))
879  res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, shape(res))
880  end select
881  type is (integer(i4_kind))
882  allocate(integer(i4_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
883  select type(res)
884  type is (integer(i4_kind))
885  res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, shape(res))
886  end select
887  end select
888 
889 end subroutine get_remapped_diurnal_data
890 
891 !> @brief Determine if there is any data to write (i.e send_data has been called)
892 !! @return .true. if there is data to write
893 function is_there_data_to_write(this) &
894  result(res)
895  class(fmsdiagoutputbuffer_type), intent(in) :: this !< Buffer object
896 
897  logical :: res
898 
899  if (allocated(this%send_data_called)) then
900  res = this%send_data_called
901  else
902  res = .false.
903  endif
904 end function
905 
906 !> @brief Determine if it is time to finish the reduction method
907 !! @return .true. if it is time to finish the reduction method
908 function is_time_to_finish_reduction(this, end_time) &
909  result(res)
910  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< Buffer object
911  type(time_type), optional, intent(in) :: end_time !< The time at the end of the run
912 
913  logical :: res
914 
915  res = .false.
916  if (this%time >= this%next_output) res = .true.
917 
918  if (present(end_time)) then
919  if (end_time >= this%next_output) res = .true.
920  endif
921 end function is_time_to_finish_reduction
922 
923 !> @brief Sets send_data_called to .true.
924 subroutine set_send_data_called(this)
925  class(fmsdiagoutputbuffer_type), intent(inout) :: this !< Buffer object
926 
927  this%send_data_called = .true.
928 end subroutine set_send_data_called
929 #endif
type(time_type) function get_base_time()
gets the module variable base_time
Definition: diag_data.F90:511
integer, parameter time_min
The reduction method is min value.
Definition: diag_data.F90:117
integer, parameter time_max
The reduction method is max value.
Definition: diag_data.F90:118
integer, parameter r8
Supported type/kind of the variable.
Definition: diag_data.F90:83
Write data to a defined field within a file Example usage:
Definition: fms2_io.F90:261
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 ...
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:406
Error handler.
Definition: mpp.F90:381
integer function, public get_ticks_per_second()
Returns the number of ticks per second.
subroutine, public get_time(Time, seconds, days, ticks, err_msg)
Returns days and seconds ( < 86400 ) corresponding to a time. err_msg should be checked for any error...
Type to represent amounts of time. Implemented as seconds and days to allow for larger intervals.
Does the time_none reduction method. See include/fms_diag_reduction_methods.inc.
Contains buffer types and routines for the diag manager.
Data structure holding a 3D bounding box. It is commonlyused to represent the interval bounds or limi...