FMS  2025.04
Flexible Modeling System
gather_data_bc.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 !> @file
19 !> @brief Routines for the @ref gather_data_bc interface
20 
21 !> @addtogroup netcdf_io_mod
22 !> @{
23 
24 !> @brief gathers the 2d vdata from all of the relevant pes into the root_pe and saves it to a
25 !! buffer.
26 subroutine gather_data_bc_2d(fileobj, vdata, bc_info)
27  class(fmsnetcdffile_t), intent(inout) :: fileobj !< Fms2io netcdf fileobj
28  class(*), dimension(:,:), intent(in) :: vdata !< Data to be gather
29  type(bc_information), intent(inout) :: bc_info !< information about the boundary condition variable
30 
31  integer :: i_glob !< Size of the global domain in x
32  integer :: j_glob !< Size of the global domain in y
33 
34  integer :: isc, iec, jsc, jec !< current PE's indices
35  integer :: i1, i2, j1, j2 !< current PE's indices relative to the global domain
36  integer :: i_add, j_add !< current PE's shift
37 
38  real(kind=r4_kind), dimension(:,:), allocatable, target :: global_buf_r4_kind !< buffer with a data gathered
39  real(kind=r4_kind), dimension(:,:), allocatable, target :: local_buf_r4_kind !< current PEs data
40  real(kind=r8_kind), dimension(:,:), allocatable, target :: global_buf_r8_kind !< buffer with a data gathered
41  real(kind=r8_kind), dimension(:,:), allocatable, target :: local_buf_r8_kind !< current PEs data
42 
43  integer(kind=i8_kind) :: chksum_val !< Checksum value calculated from the gathered data
44  character(len=32) :: chksum !< Cheksum value converted to a string
45 
46  !> Set the indices
47  isc = bc_info%indices(1)
48  iec = bc_info%indices(2)
49  jsc = bc_info%indices(3)
50  jec = bc_info%indices(4)
51 
52  !> This is the section of the PE that will actually be added to the global_buffer
53  i1 = 1 + bc_info%x_halo
54  i2 = i1 + (iec-isc)
55  j1 = 1 + bc_info%y_halo
56  j2 = j1 + (jec-jsc)
57 
58  !> Set up index shifts for global array
59  i_add = bc_info%ishift
60  j_add = bc_info%jshift
61 
62  !> Only root allocates a global_buffer that will be written
63  if (fileobj%is_root) then
64  i_glob = bc_info%global_size(1)
65  j_glob = bc_info%global_size(2)
66  endif
67 
68  !> Gather the data and calculate the checksum for the resulting array.
69  select type(vdata)
70  type is (real(kind=r4_kind))
71  !> If the fileobj's root does not have data for this variable
72  if (fileobj%is_root .and. .not. bc_info%data_on_file_root) then
73  !> Allocate global_buf_r4_kind to be one size bigger, global_buf_r4_kind(i_glob+1,,:) is just dummy data
74  allocate(global_buf_r4_kind(i_glob+1, j_glob))
75  !> Allocate a temp local buffer to the fileobj's root. This is needed because the data needs to be send
76  !! to the fileobj's root, but because fileobj's root doesn't have any data, we just create dummy data and
77  !! not write it later
78  allocate(local_buf_r4_kind(1,1))
79  local_buf_r4_kind = 0.
80  isc = 1+i_glob; i_add=0; iec=1+i_glob; jsc=j_glob; j_add=0; jec=j_glob
81  i1=1; i2=1; j1=1; j2=1
82  else
83  !! In this case there is data in fileobj's root, so there is no need for the dummy data
84  if(fileobj%is_root) allocate(global_buf_r4_kind(i_glob, j_glob))
85  allocate(local_buf_r4_kind(size(vdata,1), size(vdata,2)))
86  local_buf_r4_kind = vdata
87  endif
88 
89  call mpp_gather(isc+i_add, iec+i_add, jsc+j_add, jec+j_add, bc_info%pelist, &
90  local_buf_r4_kind(i1:i2,j1:j2), &
91  global_buf_r4_kind, fileobj%is_root)
92 
93  deallocate(local_buf_r4_kind)
94  !> If you are on fileobj's root, calculate the checksum and save the gathered data in a buffer
95  if (fileobj%is_root) then
96  chksum_val = mpp_chksum(global_buf_r4_kind(1:i_glob,1:j_glob), (/mpp_pe()/))
97  allocate(bc_info%globaldata2d_r4(i_glob, j_glob))
98  bc_info%globaldata2d_r4=global_buf_r4_kind(1:i_glob,1:j_glob)
99  deallocate(global_buf_r4_kind)
100  endif
101  type is (real(kind=r8_kind))
102  !> If the fileobj's root does not have data for this variable
103  if (fileobj%is_root .and. .not. bc_info%data_on_file_root) then
104  !> Allocate global_buf_r8_kind to be one size bigger, global_buf_r8_kind(i_glob+1,,:) is just dummy data
105  allocate(global_buf_r8_kind(i_glob+1, j_glob))
106  !> Allocate a temp local buffer to the fileobj's root. This is needed because the data needs to be send
107  !! to the fileobj's root, but because fileobj's root doesn't have any data, we just create dummy data and
108  !! not write it later
109  allocate(local_buf_r8_kind(1,1))
110  local_buf_r8_kind = 0.
111  isc = 1+i_glob; i_add=0; iec=1+i_glob; jsc=j_glob; j_add=0; jec=j_glob
112  i1=1; i2=1; j1=1; j2=1
113  else
114  !! In this case there is data in fileobj's root, so there is no need for the dummy data
115  if(fileobj%is_root) allocate(global_buf_r8_kind(i_glob, j_glob))
116  allocate(local_buf_r8_kind(size(vdata,1), size(vdata,2)))
117  local_buf_r8_kind = vdata
118  endif
119 
120  call mpp_gather(isc+i_add, iec+i_add, jsc+j_add, jec+j_add, bc_info%pelist, &
121  local_buf_r8_kind(i1:i2,j1:j2), &
122  global_buf_r8_kind, fileobj%is_root)
123  deallocate(local_buf_r8_kind)
124  !> If you are on fileobj's root, calculate the checksum and save the gathered data in a buffer
125  if (fileobj%is_root) then
126  chksum_val = mpp_chksum(global_buf_r8_kind(1:i_glob,1:j_glob), (/mpp_pe()/))
127  allocate(bc_info%globaldata2d_r8(i_glob, j_glob))
128  bc_info%globaldata2d_r8=global_buf_r8_kind(1:i_glob,1:j_glob)
129  deallocate(global_buf_r8_kind)
130  endif
131  class default
132  call error("gather_data_bc_2d: unsupported type. Currently only r8_kind and r4_kinds are supported")
133  end select
134 
135  !> Save the checksum, so you can write it later
136  if (fileobj%is_root) then
137  chksum = ""
138  write(chksum, "(Z16)") chksum_val
139  bc_info%chksum = chksum
140  endif
141 
142 end subroutine gather_data_bc_2d
143 
144 !> @brief gathers the 2d vdata from all of the relevant pes into the root_pe and saves it to a
145 !! buffer.
146 subroutine gather_data_bc_3d(fileobj, vdata, bc_info)
147  class(fmsnetcdffile_t), intent(inout) :: fileobj !< Fms2io netcdf fileobj
148  class(*), dimension(:,:,:), intent(in) :: vdata !< Data to be gather
149  type(bc_information), intent(inout) :: bc_info !< information about the boundary condition variable
150 
151  integer :: i_glob !< Size of the global domain in x
152  integer :: j_glob !< Size of the global domain in y
153  integer :: k_glob !< Size of the z axis
154 
155  integer :: isc, iec, jsc, jec !< current PE's indices
156  integer :: i1, i2, j1, j2 !< current PE's indices relative to the global domain
157  integer :: i_add, j_add !< current PE's shift
158 
159  real(kind=r4_kind), dimension(:,:,:), allocatable, target :: global_buf_r4_kind !< buffer with a data gathered
160  real(kind=r4_kind), dimension(:,:,:), allocatable, target :: local_buf_r4_kind !< current PEs data
161  real(kind=r8_kind), dimension(:,:,:), allocatable, target :: global_buf_r8_kind !< buffer with a data gathered
162  real(kind=r8_kind), dimension(:,:,:), allocatable, target :: local_buf_r8_kind !< current PEs data
163 
164  integer(kind=i8_kind) :: chksum_val !< Checksum value calculated from the gathered data
165  character(len=32) :: chksum !< Cheksum value converted to a string
166 
167  !> Set the indices
168  isc = bc_info%indices(1)
169  iec = bc_info%indices(2)
170  jsc = bc_info%indices(3)
171  jec = bc_info%indices(4)
172 
173  !> This is the section of the PE that will actually be added to the global_buffer
174  i1 = 1 + bc_info%x_halo
175  i2 = i1 + (iec-isc)
176  j1 = 1 + bc_info%y_halo
177  j2 = j1 + (jec-jsc)
178 
179  !> Set up index shifts for global array
180  i_add = bc_info%ishift
181  j_add = bc_info%jshift
182 
183  !> Allocate a global_buffer that will be written
184  if (fileobj%is_root) then
185  i_glob = bc_info%global_size(1)
186  j_glob = bc_info%global_size(2)
187  endif
188 
189  k_glob=bc_info%global_size(3)
190  !> Gather the data and calculate the checksum for the resulting array.
191  select type(vdata)
192  type is (real(kind=r4_kind))
193  !> If the fileobj's root does not have data for this variable
194  if (fileobj%is_root .and. .not. bc_info%data_on_file_root) then
195  !> Allocate global_buf_r8_kind to be one size bigger, global_buf_r8_kind(i_glob+1,,:) is just dummy data
196  allocate(global_buf_r4_kind(i_glob+1, j_glob, bc_info%global_size(3)))
197  !> Allocate a temp local buffer to the fileobj's root. This is needed because the data needs to be send
198  !! to the fileobj's root, but because fileobj's root doesn't have any data, we just create dummy data and
199  !! not write it later
200  allocate(local_buf_r4_kind(1,1,1))
201  local_buf_r4_kind = 0.
202  isc = 1+i_glob; i_add=0; iec=1+i_glob; jsc=j_glob; j_add=0; jec=j_glob
203  i1=1; i2=1; j1=1; j2=1
204  else
205  !! In this case there is data in fileobj's root, so there is no need for the dummy data
206  if(fileobj%is_root) allocate(global_buf_r4_kind(i_glob, j_glob, k_glob))
207  allocate(local_buf_r4_kind(size(vdata,1), size(vdata,2), size(vdata,3)))
208  local_buf_r4_kind = vdata
209  endif
210 
211  call mpp_gather(isc+i_add, iec+i_add, jsc+j_add, jec+j_add, k_glob, bc_info%pelist, &
212  local_buf_r4_kind(i1:i2,j1:j2,:), &
213  global_buf_r4_kind, fileobj%is_root)
214  deallocate(local_buf_r4_kind)
215  !> If you are on fileobj's root, calculate the checksum and save the gathered data in a buffer
216  if (fileobj%is_root) then
217  chksum_val = mpp_chksum(global_buf_r4_kind(1:i_glob,1:j_glob, :), (/mpp_pe()/))
218  allocate(bc_info%globaldata3d_r4(i_glob, j_glob, bc_info%global_size(3)))
219  bc_info%globaldata3d_r4=global_buf_r4_kind(1:i_glob,1:j_glob,:)
220  deallocate(global_buf_r4_kind)
221  endif
222  type is (real(kind=r8_kind))
223  !> If the fileobj's root does not have data for this variable
224  if (fileobj%is_root .and. .not. bc_info%data_on_file_root) then
225  !> Allocate global_buf_r8_kind to be one size bigger, global_buf_r8_kind(i_glob+1,,:) is just dummy data
226  allocate(global_buf_r8_kind(i_glob+1, j_glob, bc_info%global_size(3)))
227  !> Allocate a temp local buffer to the fileobj's root. This is needed because the data needs to be send
228  !! to the fileobj's root, but because fileobj's root doesn't have any data, we just create dummy data and
229  !! not write it later
230  allocate(local_buf_r8_kind(1,1,1))
231  local_buf_r8_kind = 0.
232  isc = 1+i_glob; i_add=0; iec=1+i_glob; jsc=j_glob; j_add=0; jec=j_glob
233  i1=1; i2=1; j1=1; j2=1
234  else
235  !! In this case there is data in fileobj's root, so there is no need for the dummy data
236  if(fileobj%is_root) allocate(global_buf_r8_kind(i_glob, j_glob, k_glob))
237  allocate(local_buf_r8_kind(size(vdata,1), size(vdata,2), size(vdata,3)))
238  local_buf_r8_kind = vdata
239  endif
240 
241  call mpp_gather(isc+i_add, iec+i_add, jsc+j_add, jec+j_add, k_glob, bc_info%pelist, &
242  local_buf_r8_kind(i1:i2,j1:j2,:), &
243  global_buf_r8_kind, fileobj%is_root)
244  deallocate(local_buf_r8_kind)
245  !> If you are on fileobj's root, calculate the checksum and save the gathered data in a buffer
246  if (fileobj%is_root) then
247  chksum_val = mpp_chksum(global_buf_r8_kind(1:i_glob,1:j_glob, :), (/mpp_pe()/))
248  allocate(bc_info%globaldata3d_r8(i_glob, j_glob, bc_info%global_size(3)))
249  bc_info%globaldata3d_r8=global_buf_r8_kind(1:i_glob,1:j_glob,:)
250  deallocate(global_buf_r8_kind)
251  endif
252  class default
253  call error("gather_data_bc_3d: unsupported type. Currently only r8_kind and r4_kinds are supported")
254  end select
255 
256  !> Save the checksum, so you can write it later
257  if (fileobj%is_root) then
258  chksum = ""
259  write(chksum, "(Z16)") chksum_val
260  bc_info%chksum = chksum
261  endif
262 
263 end subroutine gather_data_bc_3d
264 !> @}
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:406
subroutine gather_data_bc_3d(fileobj, vdata, bc_info)
gathers the 2d vdata from all of the relevant pes into the root_pe and saves it to a buffer.
subroutine gather_data_bc_2d(fileobj, vdata, bc_info)
gathers the 2d vdata from all of the relevant pes into the root_pe and saves it to a buffer.