FMS  2024.03
Flexible Modeling System
scatter_data_bc.inc
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 !> @file
20 !> @brief Routines for the @ref gather_data_bc interface
21 
22 !> @addtogroup netcdf_io_mod
23 !> @{
24 
25 !< Root pe reads the data in the netcdf file and scatters it to the other pes
26 subroutine scatter_data_bc_2d(fileobj, varname, vdata, bc_info, unlim_dim_level, ignore_checksum)
27  class(fmsnetcdffile_t), intent(inout) :: fileobj !< Fms2io netcdf fileobj
28  character(len=*), intent(in) :: varname !< Variable name.
29  class(*), dimension(:,:), intent(inout) :: vdata !< scattered data
30  type(bc_information), intent(inout) :: bc_info !< information about the boundary condition variable
31  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension level.
32  logical, intent(in), optional :: ignore_checksum !< Checksum data integrity flag.
33 
34  real(kind=r8_kind), dimension(:,:), allocatable :: global_buf_r8_kind !< buffer with a data read in from file
35  real(kind=r8_kind), dimension(:,:), allocatable :: local_buf_r8_kind !< current pe's portion of the data
36  real(kind=r4_kind), dimension(:,:), allocatable :: global_buf_r4_kind !< buffer with a data read in from file
37  real(kind=r4_kind), dimension(:,:), allocatable :: local_buf_r4_kind !< current pe's portion of the data
38 
39  integer(kind=i8_kind) :: chksum_val !< Checksum value calculated from the read data
40  character(len=32) :: chksum !< Cheksum value converted to a string
41  character(len=32) :: chksum_read !< String checksum read in from file
42  logical :: chksum_ignore = .false. !< local variable for data integrity checks
43  !! default: .FALSE. - checks enabled
44 
45  integer :: isc, iec, jsc, jec !< current PE's indices
46  integer :: i1, i2, j1, j2 !< current PE's indices relative to the global domain
47  integer :: i_add, j_add !< current PE's shift
48 
49  integer :: i_glob !< Size of the global domain in x
50  integer :: j_glob !< Size of the global domain in y
51 
52  if (PRESENT(ignore_checksum)) chksum_ignore = ignore_checksum
53 
54  !> Set the indices
55  isc = bc_info%indices(1)
56  iec = bc_info%indices(2)
57  jsc = bc_info%indices(3)
58  jec = bc_info%indices(4)
59 
60  !> This is the section of the PE that will actually be read from the global_buffer
61  i1 = 1 + bc_info%x_halo
62  i2 = i1 + (iec-isc)
63  j1 = 1 + bc_info%y_halo
64  j2 = j1 + (jec-jsc)
65 
66  !> Set up index shifts for global array
67  i_add = bc_info%ishift
68  j_add = bc_info%jshift
69 
70  if (fileobj%is_root) then
71  i_glob = bc_info%global_size(1)
72  j_glob = bc_info%global_size(2)
73  endif
74 
75  select type(vdata)
76  type is (real(kind=r4_kind))
77  if (fileobj%is_root) then
78  !> If you are the file root, read in the data
79  allocate(global_buf_r4_kind(i_glob, j_glob))
80 
81  call netcdf_read_data(fileobj, varname, &
82  global_buf_r4_kind, &
83  unlim_dim_level=unlim_dim_level, &
84  broadcast=.false.)
85  !> If the checksum exists read it and compare it with the calculated from the data that was read
86  if (.not.chksum_ignore .and. variable_att_exists(fileobj, varname, "checksum", broadcast=.false.)) then
87  call get_variable_attribute(fileobj, varname, "checksum", chksum_read, broadcast=.false.)
88  chksum_val = mpp_chksum(global_buf_r4_kind, (/mpp_pe()/))
89  chksum = ""
90  write(chksum, "(Z16)") chksum_val
91  if (.not. string_compare(adjustl(chksum), adjustl(chksum_read))) then
92  call error("scatter_data_bc_2d: "//trim(varname)//" chksum_in:"//(chksum)//" chksum_out:"// &
93  & (chksum_read))
94  endif
95  endif
96  endif
97 
98  !> If the fileobj's root is not the same as the variable's root, then no data will be read for the
99  !! file root
100  if (fileobj%is_root .and. .not. bc_info%data_on_file_root) then
101  allocate(local_buf_r4_kind(1,1))
102  local_buf_r4_kind = 0.
103  isc = 1; i_add=0; iec=1; jsc=1; j_add=0; jec=1
104  i1=1; i2=1; j1=1; j2=1
105  else
106  allocate(local_buf_r4_kind(i1:i2,j1:j2))
107  endif
108 
109  !> Scatter the data to the other PE
110  call mpp_scatter(isc+i_add, iec+i_add, jsc+j_add, jec+j_add, bc_info%pelist, &
111  local_buf_r4_kind(i1:i2,j1:j2), global_buf_r4_kind, fileobj%is_root)
112 
113  !> Return if fileobj's root is not the same as the variable's root
114  if (fileobj%is_root .and. .not. bc_info%data_on_file_root) then
115  deallocate(local_buf_r4_kind)
116  return
117  endif
118 
119  !> Set vdata to equal the scattered data
120  vdata(i1:i2,j1:j2) = local_buf_r4_kind(i1:i2,j1:j2)
121 
122  deallocate(local_buf_r4_kind)
123  if (fileobj%is_root) deallocate(global_buf_r4_kind)
124 
125  type is (real(kind=r8_kind))
126  if (fileobj%is_root) then
127  !> If you are the file root, read in the data
128  allocate(global_buf_r8_kind(i_glob, j_glob))
129 
130  call netcdf_read_data(fileobj, varname, &
131  global_buf_r8_kind, &
132  unlim_dim_level=unlim_dim_level, &
133  broadcast=.false.)
134  !> If the checksum exists read it and compare it with the calculated from the data that was read
135  if (.not.chksum_ignore .and. variable_att_exists(fileobj, varname, "checksum", broadcast=.false.)) then
136  call get_variable_attribute(fileobj, varname, "checksum", chksum_read, broadcast=.false.)
137  chksum_val = mpp_chksum(global_buf_r8_kind, (/mpp_pe()/))
138  chksum = ""
139  write(chksum, "(Z16)") chksum_val
140  if (.not. string_compare(adjustl(chksum), adjustl(chksum_read))) then
141  call error("scatter_data_bc_2d: "//trim(varname)//" chksum_in:"//chksum//" chksum_out:"//chksum_read)
142  endif
143  endif
144  endif
145 
146  !> If the fileobj's root is not the same as the variable's root, then no data will be read for the
147  !! file root
148  if (fileobj%is_root .and. .not. bc_info%data_on_file_root) then
149  allocate(local_buf_r8_kind(1,1))
150  local_buf_r8_kind = 0.
151  isc = 1; i_add=0; iec=1; jsc=1; j_add=0; jec=1
152  i1=1; i2=1; j1=1; j2=1
153  else
154  allocate(local_buf_r8_kind(i1:i2,j1:j2))
155  endif
156 
157  !> Scatter the data to the other PE
158  call mpp_scatter(isc+i_add, iec+i_add, jsc+j_add, jec+j_add, bc_info%pelist, &
159  local_buf_r8_kind(i1:i2,j1:j2), global_buf_r8_kind, fileobj%is_root)
160 
161  !> Return if fileobj's root is not the same as the variable's root
162  if (fileobj%is_root .and. .not. bc_info%data_on_file_root) then
163  deallocate(local_buf_r8_kind)
164  return
165  endif
166 
167  !> Set vdata to equal the scattered data
168  vdata(i1:i2,j1:j2) = local_buf_r8_kind(i1:i2,j1:j2)
169 
170  deallocate(local_buf_r8_kind)
171  if (fileobj%is_root) deallocate(global_buf_r8_kind)
172 
173  class default
174  call error("scatter_data_bc_2d: unsupported type. Currently only r8_kind and r4_kinds are supported")
175  end select
176 end subroutine scatter_data_bc_2d
177 
178 !< Root pe reads the data in the netcdf file and scatters it to the other pes
179 subroutine scatter_data_bc_3d(fileobj, varname, vdata, bc_info, unlim_dim_level, ignore_checksum)
180  class(fmsnetcdffile_t), intent(inout) :: fileobj !< Fms2io netcdf fileobj
181  character(len=*), intent(in) :: varname !< Variable name.
182  class(*), dimension(:,:,:), intent(inout) :: vdata !< scattered data
183  type(bc_information), intent(inout) :: bc_info !< information about the boundary condition variable
184  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension level.
185  logical, intent(in), optional :: ignore_checksum !< Checksum data integrity flag.
186 
187  real(kind=r8_kind), dimension(:,:,:), allocatable :: global_buf_r8_kind !< buffer with a data read in from file
188  real(kind=r8_kind), dimension(:,:,:), allocatable :: local_buf_r8_kind !< current pe's portion of the data
189  real(kind=r4_kind), dimension(:,:,:), allocatable :: global_buf_r4_kind !< buffer with a data read in from file
190  real(kind=r4_kind), dimension(:,:,:), allocatable :: local_buf_r4_kind !< current pe's portion of the data
191 
192  integer(kind=i8_kind) :: chksum_val !< Checksum value calculated from the read data
193  character(len=32) :: chksum !< Cheksum value converted to a string
194  character(len=32) :: chksum_read !< String checksum read in from file
195  logical :: chksum_ignore = .false. !< local variable for data integrity checks
196  !! default: .FALSE. - checks enabled
197 
198  integer :: isc, iec, jsc, jec !< current PE's indices
199  integer :: i1, i2, j1, j2 !< current PE's indices relative to the global domain
200  integer :: i_add, j_add !< current PE's shift
201 
202  integer :: i_glob !< Size of the global domain in x
203  integer :: j_glob !< Size of the global domain in y
204  integer :: k_glob !< Size of the z axis
205 
206  if (PRESENT(ignore_checksum)) chksum_ignore = ignore_checksum
207 
208  !> Set the indices
209  isc = bc_info%indices(1)
210  iec = bc_info%indices(2)
211  jsc = bc_info%indices(3)
212  jec = bc_info%indices(4)
213 
214  !> This is the section of the PE that will actually be read from the global_buffer
215  i1 = 1 + bc_info%x_halo
216  i2 = i1 + (iec-isc)
217  j1 = 1 + bc_info%y_halo
218  j2 = j1 + (jec-jsc)
219 
220  !> Set up index shifts for global array
221  i_add = bc_info%ishift
222  j_add = bc_info%jshift
223 
224  if (fileobj%is_root) then
225  i_glob = bc_info%global_size(1)
226  j_glob = bc_info%global_size(2)
227  endif
228 
229  k_glob = bc_info%global_size(3)
230 
231  select type(vdata)
232  type is (real(kind=r4_kind))
233  if (fileobj%is_root) then
234  !> If you are the file root, read in the data
235  allocate(global_buf_r4_kind(i_glob, j_glob, k_glob))
236 
237  call netcdf_read_data(fileobj, varname, &
238  global_buf_r4_kind, &
239  unlim_dim_level=unlim_dim_level, &
240  broadcast=.false.)
241  !> If the checksum exists read it and compare it with the calculated from the data that was read
242  if (.not.chksum_ignore .and. variable_att_exists(fileobj, varname, "checksum", broadcast=.false.)) then
243  call get_variable_attribute(fileobj, varname, "checksum", chksum_read, broadcast=.false.)
244  chksum_val = mpp_chksum(global_buf_r4_kind, (/mpp_pe()/))
245  chksum = ""
246  write(chksum, "(Z16)") chksum_val
247  if (.not. string_compare(adjustl(chksum), adjustl(chksum_read))) then
248  call error("scatter_data_bc_2d: "//trim(varname)//" chksum_in:"//(chksum)//" chksum_out:"// &
249  & (chksum_read))
250  endif
251  endif
252  endif
253 
254  !> If the fileobj's root is not the same as the variable's root, then no data will be read for the
255  !! file root
256  if (fileobj%is_root .and. .not. bc_info%data_on_file_root) then
257  allocate(local_buf_r4_kind(1,1,k_glob))
258  local_buf_r4_kind = 0.
259  isc = 1; i_add=0; iec=1; jsc=1; j_add=0; jec=1
260  i1=1; i2=1; j1=1; j2=1
261  else
262  allocate(local_buf_r4_kind(i1:i2,j1:j2,k_glob))
263  endif
264 
265  !> Scatter the data to the other PE
266  call mpp_scatter(isc+i_add, iec+i_add, jsc+j_add, jec+j_add, k_glob, bc_info%pelist, &
267  local_buf_r4_kind(i1:i2,j1:j2,:), global_buf_r4_kind, fileobj%is_root)
268 
269  !> Return if fileobj's root is not the same as the variable's root
270  if (fileobj%is_root .and. .not. bc_info%data_on_file_root) then
271  deallocate(local_buf_r4_kind)
272  return
273  endif
274 
275  !> Set vdata to equal the scattered data
276  vdata(i1:i2,j1:j2,:) = local_buf_r4_kind(i1:i2,j1:j2,:)
277 
278  deallocate(local_buf_r4_kind)
279  if (fileobj%is_root) deallocate(global_buf_r4_kind)
280 
281  type is (real(kind=r8_kind))
282  if (fileobj%is_root) then
283  !> If you are the file root, read in the data
284  allocate(global_buf_r8_kind(i_glob, j_glob, k_glob))
285 
286  call netcdf_read_data(fileobj, varname, &
287  global_buf_r8_kind, &
288  unlim_dim_level=unlim_dim_level, &
289  broadcast=.false.)
290  !> If the checksum exists read it and compare it with the calculated from the data that was read
291  if (.not.chksum_ignore .and. variable_att_exists(fileobj, varname, "checksum", broadcast=.false.)) then
292  call get_variable_attribute(fileobj, varname, "checksum", chksum_read, broadcast=.false.)
293  chksum_val = mpp_chksum(global_buf_r8_kind, (/mpp_pe()/))
294  chksum = ""
295  write(chksum, "(Z16)") chksum_val
296  if (.not. string_compare(adjustl(chksum), adjustl(chksum_read))) then
297  call error("scatter_data_bc_2d: "//trim(varname)//" chksum_in:"//chksum//" chksum_out:"//chksum_read)
298  endif
299  endif
300  endif
301 
302  !> If the fileobj's root is not the same as the variable's root, then no data will be read for the
303  !! file root
304  if (fileobj%is_root .and. .not. bc_info%data_on_file_root) then
305  allocate(local_buf_r8_kind(1,1,k_glob))
306  local_buf_r8_kind = 0.
307  isc = 1; i_add=0; iec=1; jsc=1; j_add=0; jec=1
308  i1=1; i2=1; j1=1; j2=1;
309  else
310  allocate(local_buf_r8_kind(i1:i2,j1:j2,k_glob))
311  endif
312 
313  !> Scatter the data to the other PE
314  call mpp_scatter(isc+i_add, iec+i_add, jsc+j_add, jec+j_add, k_glob, bc_info%pelist, &
315  local_buf_r8_kind(i1:i2,j1:j2,:), global_buf_r8_kind, fileobj%is_root)
316 
317  !> Return if fileobj's root is not the same as the variable's root
318  if (fileobj%is_root .and. .not. bc_info%data_on_file_root) then
319  deallocate(local_buf_r8_kind)
320  return
321  endif
322 
323  !> Set vdata to equal the scattered data
324  vdata(i1:i2,j1:j2,:) = local_buf_r8_kind(i1:i2,j1:j2,:)
325 
326  deallocate(local_buf_r8_kind)
327  if (fileobj%is_root) deallocate(global_buf_r8_kind)
328 
329  class default
330  call error("scatter_data_bc_3d: unsupported type. Currently only r8_kind and r4_kinds are supported")
331  end select
332 end subroutine scatter_data_bc_3d
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407
subroutine scatter_data_bc_2d(fileobj, varname, vdata, bc_info, unlim_dim_level, ignore_checksum)
subroutine scatter_data_bc_3d(fileobj, varname, vdata, bc_info, unlim_dim_level, ignore_checksum)