FMS  2025.04
Flexible Modeling System
fms_diag_input_buffer.inc
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 
19 !> @brief Appends the input_data_buffer and the mask (only when the mask is set to .True.)
20 subroutine append_data_buffer_(mask_out, mask_in, data_out, data_in)
21  logical, intent(inout) :: mask_out(:,:,:,:) !< Mask currently in the input_data_buffer
22  logical, intent(in) :: mask_in(:,:,:,:) !< Mask passed in to send_data
23  real(FMS_TRM_KIND_), intent(inout) :: data_out(:,:,:,:) !< Data currently in the input_data_buffer
24  real(FMS_TRM_KIND_), intent(in) :: data_in(:,:,:,:) !< Data passed in to send_data
25 
26  integer :: i, j, k, l !< For looping through the input_data_buffer
27 
28  do l = 1, size(data_out, 4)
29  do k = 1, size(data_out, 3)
30  do j = 1, size(data_out, 2)
31  do i = 1, size(data_out, 1)
32  if (mask_in(i,j,k,l)) then
33  mask_out(i,j,k,l) = .true.
34  data_out(i,j,k,l) = data_in(i,j,k,l)
35  endif
36  enddo
37  enddo
38  enddo
39  enddo
40 
41 end subroutine append_data_buffer_
42 
43 !> @brief Sums the data in the input_data_buffer
44 subroutine sum_data_buffer_(mask, data_out, data_in, counter, var_is_masked)
45  logical, intent(in) :: mask(:,:,:,:) !< Mask passed into send_data
46  real(FMS_TRM_KIND_), intent(inout) :: data_out(:,:,:,:) !< Data currently saved in the input_data_buffer
47  real(FMS_TRM_KIND_), intent(in) :: data_in(:,:,:,:) !< Data passed into send_data
48  integer, intent(inout) :: counter(:,:,:,:) !< Number of times data has been summed
49  logical, intent(in) :: var_is_masked !< .True. if the variable is masked
50 
51  if (var_is_masked) then
52  where (mask)
53  data_out = data_out + data_in
54  endwhere
55  else
56  data_out = data_out + data_in
57  endif
58 
59  counter = counter + 1
60 end subroutine SUM_DATA_BUFFER_