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