30 use constants_mod,
only: seconds_per_day
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
38 use fms_diag_time_utils_mod,
only: diag_time_inc
47 integer(i4_kind) :: buffer_type
48 class(*),
allocatable :: buffer(:,:,:,:,:)
50 integer,
allocatable :: buffer_dims(:)
51 real(r8_kind),
allocatable :: weight_sum(:,:,:,:)
54 integer,
allocatable :: num_elements(:)
55 integer,
allocatable :: axis_ids(:)
58 logical :: done_with_math
59 integer :: diurnal_sample_size = -1
61 integer :: diurnal_section= -1
63 logical,
allocatable :: send_data_called
64 integer :: unlmited_dimension
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
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
113 public :: fms_diag_output_buffer_init
121 logical function fms_diag_output_buffer_init(buffobjs, buff_list_size)
124 integer,
intent(in) :: buff_list_size
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
134 subroutine set_buffer_id(this, id)
136 integer,
intent(in) :: id
139 end subroutine set_buffer_id
142 subroutine flush_buffer(this)
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
158 subroutine allocate_buffer(this, buff_type, ndim, buff_sizes, mask_variant, field_name, diurnal_samples)
160 class(*),
intent(in) :: buff_type
161 integer,
intent(in) :: ndim
162 integer,
intent(in) :: buff_sizes(4)
163 logical,
intent(in) :: mask_variant
164 character(len=*),
intent(in) :: field_name
165 integer,
intent(in) :: diurnal_samples
169 n_samples = max(1, diurnal_samples)
170 call this%set_diurnal_sample_size(n_samples)
173 if(
allocated(this%buffer))
call mpp_error(fatal,
"allocate_buffer: buffer already allocated for field:" // &
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), &
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), &
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), &
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), &
191 this%buffer_type =
r8
194 "The buff_type value passed to allocate a buffer is not a r8, r4, i8, or i4" // &
195 "for field:" // field_name, fatal)
197 if (mask_variant)
then
198 allocate(this%weight_sum(buff_sizes(1),buff_sizes(2),buff_sizes(3),buff_sizes(4)))
200 allocate(this%weight_sum(1,1,1,1))
202 this%weight_sum = 0.0_r8_kind
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
218 subroutine get_buffer (this, buff_out, field_name)
220 class(*),
allocatable,
intent(out) :: buff_out(:,:,:,:,:)
222 character(len=*),
intent(in) :: field_name
224 integer(i4_kind) :: buff_size(5)
226 if(.not.
allocated(this%buffer))
call mpp_error(fatal,
'get_buffer: buffer not yet allocated for field:' &
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)
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)))
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)))
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)))
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)))
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)
254 subroutine initialize_buffer (this, reduction_method, field_name)
256 integer,
intent(in) :: reduction_method
257 character(len=*),
intent(in) :: field_name
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.')
262 select type(buff => this%buffer)
263 type is(real(r8_kind))
264 select case (reduction_method)
266 buff = real(min_value, kind=r8_kind)
268 buff = real(max_value, kind=r8_kind)
270 buff = real(empty, kind=r8_kind)
272 type is(real(r4_kind))
273 select case (reduction_method)
275 buff = real(min_value, kind=r4_kind)
277 buff = real(max_value, kind=r4_kind)
279 buff = real(empty, kind=r4_kind)
281 type is(
integer(i8_kind))
282 select case (reduction_method)
284 buff = int(min_value, kind=i8_kind)
286 buff = int(max_value, kind=i8_kind)
288 buff = int(empty, kind=i8_kind)
290 type is(
integer(i4_kind))
291 select case (reduction_method)
293 buff = int(min_value, kind=i4_kind)
295 buff = int(max_value, kind=i4_kind)
297 buff = int(empty, kind=i4_kind)
300 call mpp_error(fatal,
'initialize buffer_5d: buffer allocated to invalid data type, this shouldnt happen')
303 end subroutine initialize_buffer
307 function get_buffer_name(this) &
311 character(len=:),
allocatable :: rslt
313 rslt = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
314 end function get_buffer_name
317 subroutine add_axis_ids(this, axis_ids)
319 integer,
intent(in) :: axis_ids(:)
321 this%axis_ids = axis_ids
326 subroutine get_axis_ids(this, res)
328 integer,
pointer,
intent(out) :: res(:)
330 if (
allocated(this%axis_ids))
then
340 function get_field_id(this) &
346 end function get_field_id
349 subroutine set_field_id(this, field_id)
351 integer,
intent(in) :: field_id
353 this%field_id = field_id
354 end subroutine set_field_id
357 subroutine set_yaml_id(this, yaml_id)
359 integer,
intent(in) :: yaml_id
361 this%yaml_id = yaml_id
362 end subroutine set_yaml_id
365 subroutine init_buffer_time(this, time)
367 type(
time_type),
optional,
intent(in) :: time
369 if (
present(time))
then
371 this%next_output = time
374 this%next_output = this%time
376 end subroutine init_buffer_time
379 subroutine set_next_output(this, next_output, next_next_output, is_static)
381 type(
time_type),
intent(in) :: next_output
382 type(
time_type),
intent(in) :: next_next_output
383 logical,
optional,
intent(in) :: is_static
385 if (
present(is_static))
then
389 this%next_output = this%time
398 if (next_output > this%next_output)
then
399 this%next_output = next_output
401 this%next_output = next_next_output
403 end subroutine set_next_output
406 subroutine update_buffer_time(this, time)
410 if (time > this%time)
then
413 end subroutine update_buffer_time
417 function get_buffer_time(this) &
423 end function get_buffer_time
427 function is_done_with_math(this) &
432 res = this%done_with_math
433 end function is_done_with_math
436 subroutine set_done_with_math(this)
440 this%done_with_math = .true.
441 end subroutine set_done_with_math
445 function get_yaml_id(this) &
452 end function get_yaml_id
456 function get_unlim_dim(this) &
461 res = this%unlmited_dimension
462 end function get_unlim_dim
465 subroutine increase_unlim_dim(this)
468 this%unlmited_dimension = this%unlmited_dimension + 1
469 end subroutine increase_unlim_dim
472 subroutine init_buffer_unlim_dim(this)
475 this%unlmited_dimension = 0
476 end subroutine init_buffer_unlim_dim
479 subroutine write_buffer(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
481 class(fmsnetcdffile_t),
intent(in) :: fms2io_fileobj
482 integer,
optional,
intent(in) :: unlim_dim_level
483 logical,
optional,
intent(in) :: is_diurnal
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)
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.")
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())
501 end subroutine write_buffer
504 subroutine write_buffer_wrapper_netcdf(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
506 type(fmsnetcdffile_t),
intent(in) :: fms2io_fileobj
507 integer,
optional,
intent(in) :: unlim_dim_level
508 logical,
optional,
intent(in) :: is_diurnal
510 character(len=:),
allocatable :: varname
511 logical :: using_diurnal
512 class(*),
allocatable :: buff_ptr(:,:,:,:,:)
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)
519 buff_ptr = this%buffer
522 varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
523 select case(this%ndim)
525 call write_data(fms2io_fileobj, varname, buff_ptr(1,1,1,1,1), unlim_dim_level=unlim_dim_level)
527 call write_data(fms2io_fileobj, varname, buff_ptr(:,1,1,1,1), unlim_dim_level=unlim_dim_level)
529 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,1,1,1), unlim_dim_level=unlim_dim_level)
531 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,1,1), unlim_dim_level=unlim_dim_level)
533 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,1), unlim_dim_level=unlim_dim_level)
535 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,:), unlim_dim_level=unlim_dim_level)
537 end subroutine write_buffer_wrapper_netcdf
540 subroutine write_buffer_wrapper_domain(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
542 type(fmsnetcdfdomainfile_t),
intent(in) :: fms2io_fileobj
543 integer,
optional,
intent(in) :: unlim_dim_level
544 logical,
optional,
intent(in) :: is_diurnal
547 character(len=:),
allocatable :: varname
548 logical :: using_diurnal
549 class(*),
allocatable :: buff_ptr(:,:,:,:,:)
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)
556 buff_ptr = this%buffer
559 varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
560 select case(this%ndim)
562 call write_data(fms2io_fileobj, varname, buff_ptr(1,1,1,1,1), unlim_dim_level=unlim_dim_level)
564 call write_data(fms2io_fileobj, varname, buff_ptr(:,1,1,1,1), unlim_dim_level=unlim_dim_level)
566 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,1,1,1), unlim_dim_level=unlim_dim_level)
568 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,1,1), unlim_dim_level=unlim_dim_level)
570 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,1), unlim_dim_level=unlim_dim_level)
572 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,:), unlim_dim_level=unlim_dim_level)
574 end subroutine write_buffer_wrapper_domain
577 subroutine write_buffer_wrapper_u(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
579 type(fmsnetcdfunstructureddomainfile_t),
intent(in) :: fms2io_fileobj
580 integer,
optional,
intent(in) :: unlim_dim_level
581 logical,
optional,
intent(in) :: is_diurnal
584 character(len=:),
allocatable :: varname
585 logical :: using_diurnal
586 class(*),
allocatable :: buff_ptr(:,:,:,:,:)
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)
593 buff_ptr = this%buffer
596 varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
597 select case(this%ndim)
599 call write_data(fms2io_fileobj, varname, buff_ptr(1,1,1,1,1), unlim_dim_level=unlim_dim_level)
601 call write_data(fms2io_fileobj, varname, buff_ptr(:,1,1,1,1), unlim_dim_level=unlim_dim_level)
603 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,1,1,1), unlim_dim_level=unlim_dim_level)
605 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,1,1), unlim_dim_level=unlim_dim_level)
607 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,1), unlim_dim_level=unlim_dim_level)
609 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,:), unlim_dim_level=unlim_dim_level)
611 end subroutine write_buffer_wrapper_u
615 function do_time_none_wrapper(this, field_data, mask, is_masked, bounds_in, bounds_out, missing_value) &
618 class(*),
intent(in) :: field_data(:,:,:,:)
621 logical,
intent(in) :: mask(:,:,:,:)
622 logical,
intent(in) :: is_masked
623 real(kind=r8_kind),
intent(in) :: missing_value
624 character(len=50) :: 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)
634 err_msg=
"do_time_none_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
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))
642 err_msg=
"do_time_none_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
645 end function do_time_none_wrapper
649 function do_time_min_wrapper(this, field_data, mask, is_masked, bounds_in, bounds_out, missing_value) &
652 class(*),
intent(in) :: field_data(:,:,:,:)
655 logical,
intent(in) :: mask(:,:,:,:)
656 logical,
intent(in) :: is_masked
657 real(kind=r8_kind),
intent(in) :: missing_value
658 character(len=50) :: 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)
668 err_msg=
"do_time_min_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
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))
676 err_msg=
"do_time_min_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
679 end function do_time_min_wrapper
683 function do_time_max_wrapper(this, field_data, mask, is_masked, bounds_in, bounds_out, missing_value) &
686 class(*),
intent(in) :: field_data(:,:,:,:)
689 logical,
intent(in) :: mask(:,:,:,:)
690 logical,
intent(in) :: is_masked
691 real(kind=r8_kind),
intent(in) :: missing_value
692 character(len=50) :: 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)
702 err_msg=
"do_time_max_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
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))
710 err_msg=
"do_time_max_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
713 end function do_time_max_wrapper
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) &
721 class(*),
intent(in) :: field_data(:,:,:,:)
724 logical,
intent(in) :: mask(:,:,:,:)
725 logical,
intent(in) :: is_masked
726 logical,
intent(in) :: mask_variant
727 real(kind=r8_kind),
intent(in) :: missing_value
728 logical,
intent(in) :: has_missing_value
730 integer,
optional,
intent(in) :: pow_value
733 real(kind=r8_kind),
optional,
intent(in) :: weight
734 character(len=150) :: 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!"
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)
750 err_msg=
"do_time_sum_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
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!"
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)
763 err_msg=
"do_time_sum_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
766 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 function do_time_sum_wrapper
773 function diag_reduction_done_wrapper(this, reduction_method, missing_value, has_mask, mask_variant) &
776 integer,
intent(in) :: reduction_method
777 real(kind=r8_kind),
intent(in) :: missing_value
778 logical,
intent(in) :: has_mask
779 logical,
intent(in) :: mask_variant
780 character(len=51) :: err_msg
782 if(.not.
allocated(this%buffer))
return
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)
793 this%weight_sum = 0.0_r8_kind
798 pure function get_buffer_dims(this)
800 integer :: get_buffer_dims(4)
801 get_buffer_dims = this%buffer_dims(1:4)
805 pure integer function get_diurnal_sample_size(this)
807 get_diurnal_sample_size = this%diurnal_sample_size
808 end function get_diurnal_sample_size
811 subroutine set_diurnal_sample_size(this, sample_size)
813 integer,
intent(in) :: sample_size
815 this%diurnal_sample_size = sample_size
816 end subroutine set_diurnal_sample_size
820 subroutine set_diurnal_section_index(this, time)
823 integer :: seconds, days, ticks
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")
828 call get_time(time,seconds,days,ticks)
831 & * this%diurnal_sample_size/seconds_per_day) + 1
832 end subroutine set_diurnal_section_index
836 subroutine get_remapped_diurnal_data(this, res)
838 class(*),
intent(out),
allocatable :: res(:,:,:,:,:)
840 integer :: ie, je, ke, ze, de
841 integer(i4_kind) :: buff_size(5)
844 last_dim = this%ndim - 1
846 ke = 1; ze = 1; de = 1
847 select case(last_dim)
849 ie = this%buffer_dims(1); je = this%buffer_dims(5)
851 ie = this%buffer_dims(1); je = this%buffer_dims(2)
852 ke = this%buffer_dims(5)
854 ie = this%buffer_dims(1); je = this%buffer_dims(2)
855 ke = this%buffer_dims(3); ze = this%buffer_dims(5)
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))
866 type is (real(r8_kind))
867 res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, shape(res))
869 type is (real(r4_kind))
870 allocate(real(r4_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
872 type is (real(r4_kind))
873 res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, shape(res))
875 type is (
integer(i8_kind))
876 allocate(
integer(i8_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
878 type is (
integer(i8_kind))
879 res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, shape(res))
881 type is (
integer(i4_kind))
882 allocate(
integer(i4_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
884 type is (
integer(i4_kind))
885 res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, shape(res))
889 end subroutine get_remapped_diurnal_data
893 function is_there_data_to_write(this) &
899 if (
allocated(this%send_data_called))
then
900 res = this%send_data_called
908 function is_time_to_finish_reduction(this, end_time) &
911 type(
time_type),
optional,
intent(in) :: end_time
916 if (this%time >= this%next_output) res = .true.
918 if (
present(end_time))
then
919 if (end_time >= this%next_output) res = .true.
921 end function is_time_to_finish_reduction
924 subroutine set_send_data_called(this)
927 this%send_data_called = .true.
928 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