31 use constants_mod,
only: seconds_per_day
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
39 use fms_diag_time_utils_mod,
only: diag_time_inc
48 integer(i4_kind) :: buffer_type
49 class(*),
allocatable :: buffer(:,:,:,:,:)
51 integer,
allocatable :: buffer_dims(:)
52 real(r8_kind),
allocatable :: weight_sum(:,:,:,:)
55 integer,
allocatable :: num_elements(:)
56 integer,
allocatable :: axis_ids(:)
59 logical :: done_with_math
60 integer :: diurnal_sample_size = -1
62 integer :: diurnal_section= -1
64 logical,
allocatable :: send_data_called
65 integer :: unlmited_dimension
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
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
114 public :: fms_diag_output_buffer_init
122 logical function fms_diag_output_buffer_init(buffobjs, buff_list_size)
125 integer,
intent(in) :: buff_list_size
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
135 subroutine set_buffer_id(this, id)
137 integer,
intent(in) :: id
140 end subroutine set_buffer_id
143 subroutine flush_buffer(this)
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
159 subroutine allocate_buffer(this, buff_type, ndim, buff_sizes, mask_variant, field_name, diurnal_samples)
161 class(*),
intent(in) :: buff_type
162 integer,
intent(in) :: ndim
163 integer,
intent(in) :: buff_sizes(4)
164 logical,
intent(in) :: mask_variant
165 character(len=*),
intent(in) :: field_name
166 integer,
intent(in) :: diurnal_samples
170 n_samples = max(1, diurnal_samples)
171 call this%set_diurnal_sample_size(n_samples)
174 if(
allocated(this%buffer))
call mpp_error(fatal,
"allocate_buffer: buffer already allocated for field:" // &
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), &
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), &
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), &
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), &
192 this%buffer_type =
r8
195 "The buff_type value passed to allocate a buffer is not a r8, r4, i8, or i4" // &
196 "for field:" // field_name, fatal)
198 if (mask_variant)
then
199 allocate(this%weight_sum(buff_sizes(1),buff_sizes(2),buff_sizes(3),buff_sizes(4)))
201 allocate(this%weight_sum(1,1,1,1))
203 this%weight_sum = 0.0_r8_kind
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
219 subroutine get_buffer (this, buff_out, field_name)
221 class(*),
allocatable,
intent(out) :: buff_out(:,:,:,:,:)
223 character(len=*),
intent(in) :: field_name
225 integer(i4_kind) :: buff_size(5)
227 if(.not.
allocated(this%buffer))
call mpp_error(fatal,
'get_buffer: buffer not yet allocated for field:' &
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)
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)))
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)))
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)))
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)))
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)
255 subroutine initialize_buffer (this, reduction_method, field_name)
257 integer,
intent(in) :: reduction_method
258 character(len=*),
intent(in) :: field_name
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.')
263 select type(buff => this%buffer)
264 type is(real(r8_kind))
265 select case (reduction_method)
267 buff = real(min_value, kind=r8_kind)
269 buff = real(max_value, kind=r8_kind)
271 buff = real(empty, kind=r8_kind)
273 type is(real(r4_kind))
274 select case (reduction_method)
276 buff = real(min_value, kind=r4_kind)
278 buff = real(max_value, kind=r4_kind)
280 buff = real(empty, kind=r4_kind)
282 type is(
integer(i8_kind))
283 select case (reduction_method)
285 buff = int(min_value, kind=i8_kind)
287 buff = int(max_value, kind=i8_kind)
289 buff = int(empty, kind=i8_kind)
291 type is(
integer(i4_kind))
292 select case (reduction_method)
294 buff = int(min_value, kind=i4_kind)
296 buff = int(max_value, kind=i4_kind)
298 buff = int(empty, kind=i4_kind)
301 call mpp_error(fatal,
'initialize buffer_5d: buffer allocated to invalid data type, this shouldnt happen')
304 end subroutine initialize_buffer
308 function get_buffer_name(this) &
312 character(len=:),
allocatable :: rslt
314 rslt = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
315 end function get_buffer_name
318 subroutine add_axis_ids(this, axis_ids)
320 integer,
intent(in) :: axis_ids(:)
322 this%axis_ids = axis_ids
327 subroutine get_axis_ids(this, res)
329 integer,
pointer,
intent(out) :: res(:)
331 if (
allocated(this%axis_ids))
then
341 function get_field_id(this) &
347 end function get_field_id
350 subroutine set_field_id(this, field_id)
352 integer,
intent(in) :: field_id
354 this%field_id = field_id
355 end subroutine set_field_id
358 subroutine set_yaml_id(this, yaml_id)
360 integer,
intent(in) :: yaml_id
362 this%yaml_id = yaml_id
363 end subroutine set_yaml_id
366 subroutine init_buffer_time(this, time)
368 type(
time_type),
optional,
intent(in) :: time
370 if (
present(time))
then
372 this%next_output = time
375 this%next_output = this%time
377 end subroutine init_buffer_time
380 subroutine set_next_output(this, next_output, next_next_output, is_static)
382 type(
time_type),
intent(in) :: next_output
383 type(
time_type),
intent(in) :: next_next_output
384 logical,
optional,
intent(in) :: is_static
386 if (
present(is_static))
then
390 this%next_output = this%time
399 if (next_output > this%next_output)
then
400 this%next_output = next_output
402 this%next_output = next_next_output
404 end subroutine set_next_output
407 subroutine update_buffer_time(this, time)
411 if (time > this%time)
then
414 end subroutine update_buffer_time
418 function get_buffer_time(this) &
424 end function get_buffer_time
428 function is_done_with_math(this) &
433 res = this%done_with_math
434 end function is_done_with_math
437 subroutine set_done_with_math(this)
441 this%done_with_math = .true.
442 end subroutine set_done_with_math
446 function get_yaml_id(this) &
453 end function get_yaml_id
457 function get_unlim_dim(this) &
462 res = this%unlmited_dimension
463 end function get_unlim_dim
466 subroutine increase_unlim_dim(this)
469 this%unlmited_dimension = this%unlmited_dimension + 1
470 end subroutine increase_unlim_dim
473 subroutine init_buffer_unlim_dim(this)
476 this%unlmited_dimension = 0
477 end subroutine init_buffer_unlim_dim
480 subroutine write_buffer(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
482 class(fmsnetcdffile_t),
intent(in) :: fms2io_fileobj
483 integer,
optional,
intent(in) :: unlim_dim_level
484 logical,
optional,
intent(in) :: is_diurnal
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)
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.")
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())
502 end subroutine write_buffer
505 subroutine write_buffer_wrapper_netcdf(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
507 type(fmsnetcdffile_t),
intent(in) :: fms2io_fileobj
508 integer,
optional,
intent(in) :: unlim_dim_level
509 logical,
optional,
intent(in) :: is_diurnal
511 character(len=:),
allocatable :: varname
512 logical :: using_diurnal
513 class(*),
allocatable :: buff_ptr(:,:,:,:,:)
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)
520 buff_ptr = this%buffer
523 varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
524 select case(this%ndim)
526 call write_data(fms2io_fileobj, varname, buff_ptr(1,1,1,1,1), unlim_dim_level=unlim_dim_level)
528 call write_data(fms2io_fileobj, varname, buff_ptr(:,1,1,1,1), unlim_dim_level=unlim_dim_level)
530 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,1,1,1), unlim_dim_level=unlim_dim_level)
532 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,1,1), unlim_dim_level=unlim_dim_level)
534 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,1), unlim_dim_level=unlim_dim_level)
536 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,:), unlim_dim_level=unlim_dim_level)
538 end subroutine write_buffer_wrapper_netcdf
541 subroutine write_buffer_wrapper_domain(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
543 type(fmsnetcdfdomainfile_t),
intent(in) :: fms2io_fileobj
544 integer,
optional,
intent(in) :: unlim_dim_level
545 logical,
optional,
intent(in) :: is_diurnal
548 character(len=:),
allocatable :: varname
549 logical :: using_diurnal
550 class(*),
allocatable :: buff_ptr(:,:,:,:,:)
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)
557 buff_ptr = this%buffer
560 varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
561 select case(this%ndim)
563 call write_data(fms2io_fileobj, varname, buff_ptr(1,1,1,1,1), unlim_dim_level=unlim_dim_level)
565 call write_data(fms2io_fileobj, varname, buff_ptr(:,1,1,1,1), unlim_dim_level=unlim_dim_level)
567 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,1,1,1), unlim_dim_level=unlim_dim_level)
569 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,1,1), unlim_dim_level=unlim_dim_level)
571 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,1), unlim_dim_level=unlim_dim_level)
573 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,:), unlim_dim_level=unlim_dim_level)
575 end subroutine write_buffer_wrapper_domain
578 subroutine write_buffer_wrapper_u(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
580 type(fmsnetcdfunstructureddomainfile_t),
intent(in) :: fms2io_fileobj
581 integer,
optional,
intent(in) :: unlim_dim_level
582 logical,
optional,
intent(in) :: is_diurnal
585 character(len=:),
allocatable :: varname
586 logical :: using_diurnal
587 class(*),
allocatable :: buff_ptr(:,:,:,:,:)
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)
594 buff_ptr = this%buffer
597 varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
598 select case(this%ndim)
600 call write_data(fms2io_fileobj, varname, buff_ptr(1,1,1,1,1), unlim_dim_level=unlim_dim_level)
602 call write_data(fms2io_fileobj, varname, buff_ptr(:,1,1,1,1), unlim_dim_level=unlim_dim_level)
604 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,1,1,1), unlim_dim_level=unlim_dim_level)
606 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,1,1), unlim_dim_level=unlim_dim_level)
608 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,1), unlim_dim_level=unlim_dim_level)
610 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,:), unlim_dim_level=unlim_dim_level)
612 end subroutine write_buffer_wrapper_u
616 function do_time_none_wrapper(this, field_data, mask, is_masked, bounds_in, bounds_out, missing_value) &
619 class(*),
intent(in) :: field_data(:,:,:,:)
622 logical,
intent(in) :: mask(:,:,:,:)
623 logical,
intent(in) :: is_masked
624 real(kind=r8_kind),
intent(in) :: missing_value
625 character(len=50) :: 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)
635 err_msg=
"do_time_none_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
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))
643 err_msg=
"do_time_none_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
646 end function do_time_none_wrapper
650 function do_time_min_wrapper(this, field_data, mask, is_masked, bounds_in, bounds_out, missing_value) &
653 class(*),
intent(in) :: field_data(:,:,:,:)
656 logical,
intent(in) :: mask(:,:,:,:)
657 logical,
intent(in) :: is_masked
658 real(kind=r8_kind),
intent(in) :: missing_value
659 character(len=50) :: 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)
669 err_msg=
"do_time_min_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
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))
677 err_msg=
"do_time_min_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
680 end function do_time_min_wrapper
684 function do_time_max_wrapper(this, field_data, mask, is_masked, bounds_in, bounds_out, missing_value) &
687 class(*),
intent(in) :: field_data(:,:,:,:)
690 logical,
intent(in) :: mask(:,:,:,:)
691 logical,
intent(in) :: is_masked
692 real(kind=r8_kind),
intent(in) :: missing_value
693 character(len=50) :: 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)
703 err_msg=
"do_time_max_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
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))
711 err_msg=
"do_time_max_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
714 end function do_time_max_wrapper
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) &
722 class(*),
intent(in) :: field_data(:,:,:,:)
725 logical,
intent(in) :: mask(:,:,:,:)
726 logical,
intent(in) :: is_masked
727 logical,
intent(in) :: mask_variant
728 real(kind=r8_kind),
intent(in) :: missing_value
729 logical,
intent(in) :: has_missing_value
731 integer,
optional,
intent(in) :: pow_value
734 real(kind=r8_kind),
optional,
intent(in) :: weight
735 character(len=150) :: 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!"
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)
751 err_msg=
"do_time_sum_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
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!"
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)
764 err_msg=
"do_time_sum_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
767 err_msg=
"do_time_sum_wrapper::the output buffer is not a valid type, must be real(r8_kind) or real(r4_kind)"
769 end function do_time_sum_wrapper
774 function diag_reduction_done_wrapper(this, reduction_method, missing_value, has_mask, mask_variant) &
777 integer,
intent(in) :: reduction_method
778 real(kind=r8_kind),
intent(in) :: missing_value
779 logical,
intent(in) :: has_mask
780 logical,
intent(in) :: mask_variant
781 character(len=51) :: err_msg
783 if(.not.
allocated(this%buffer))
return
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)
794 this%weight_sum = 0.0_r8_kind
799 pure function get_buffer_dims(this)
801 integer :: get_buffer_dims(4)
802 get_buffer_dims = this%buffer_dims(1:4)
806 pure integer function get_diurnal_sample_size(this)
808 get_diurnal_sample_size = this%diurnal_sample_size
809 end function get_diurnal_sample_size
812 subroutine set_diurnal_sample_size(this, sample_size)
814 integer,
intent(in) :: sample_size
816 this%diurnal_sample_size = sample_size
817 end subroutine set_diurnal_sample_size
821 subroutine set_diurnal_section_index(this, time)
824 integer :: seconds, days, ticks
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")
829 call get_time(time,seconds,days,ticks)
832 & * this%diurnal_sample_size/seconds_per_day) + 1
833 end subroutine set_diurnal_section_index
837 subroutine get_remapped_diurnal_data(this, res)
839 class(*),
intent(out),
allocatable :: res(:,:,:,:,:)
841 integer :: ie, je, ke, ze, de
842 integer(i4_kind) :: buff_size(5)
845 last_dim = this%ndim - 1
847 ke = 1; ze = 1; de = 1
848 select case(last_dim)
850 ie = this%buffer_dims(1); je = this%buffer_dims(5)
852 ie = this%buffer_dims(1); je = this%buffer_dims(2)
853 ke = this%buffer_dims(5)
855 ie = this%buffer_dims(1); je = this%buffer_dims(2)
856 ke = this%buffer_dims(3); ze = this%buffer_dims(5)
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))
867 type is (real(r8_kind))
868 res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, shape(res))
870 type is (real(r4_kind))
871 allocate(real(r4_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
873 type is (real(r4_kind))
874 res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, shape(res))
876 type is (
integer(i8_kind))
877 allocate(
integer(i8_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
879 type is (
integer(i8_kind))
880 res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, shape(res))
882 type is (
integer(i4_kind))
883 allocate(
integer(i4_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
885 type is (
integer(i4_kind))
886 res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, shape(res))
890 end subroutine get_remapped_diurnal_data
894 function is_there_data_to_write(this) &
900 if (
allocated(this%send_data_called))
then
901 res = this%send_data_called
909 function is_time_to_finish_reduction(this, end_time) &
912 type(
time_type),
optional,
intent(in) :: end_time
917 if (this%time >= this%next_output) res = .true.
919 if (
present(end_time))
then
920 if (end_time >= this%next_output) res = .true.
922 end function is_time_to_finish_reduction
925 subroutine set_send_data_called(this)
928 this%send_data_called = .true.
929 end subroutine set_send_data_called
type(time_type) function get_base_time()
gets the module variable base_time
integer, parameter time_min
The reduction method is min value.
integer, parameter time_max
The reduction method is max value.
integer, parameter r8
Supported type/kind of the variable.
Write data to a defined field within a file Example usage:
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.
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...
holds an allocated buffer0-5d object