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