FMS  2025.04
Flexible Modeling System
domain_read.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 Reads domain decomposed variables and scatters data to pe's for the @ref read_data interface
20 
21 !> @addtogroup fms2_io_mod
22 !> @{
23 
24 !> @brief I/O domain root reads in a domain decomposed variable at a
25 !! specific unlimited dimension level and scatters the data to the
26 !! rest of the ranks using its I/O compute domain indices. This
27 !! routine may only be used with variables that are "domain
28 !! decomposed".
29 subroutine domain_read_0d(fileobj, variable_name, vdata, unlim_dim_level, corner)
30 
31  type(fmsnetcdfdomainfile_t), intent(in) :: fileobj !< File object.
32  character(len=*), intent(in) :: variable_name !< Variable name.
33  class(*), intent(inout) :: vdata !< Data that will
34  !! be read out
35  !! to the netcdf file.
36  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
37  !! dimension.
38  integer, intent(in), optional :: corner !< Array of starting
39  !! indices describing
40  !! where the data
41  !! will be read to.
42 
43  call netcdf_read_data(fileobj, variable_name, vdata, &
44  unlim_dim_level=unlim_dim_level, corner=corner, &
45  broadcast=.true.)
46 
47 end subroutine domain_read_0d
48 
49 
50 !> @brief I/O domain root reads in a domain decomposed variable at a
51 !! specific unlimited dimension level and scatters the data to the
52 !! rest of the ranks using its I/O compute domain indices. This
53 !! routine may only be used with variables that are "domain
54 !! decomposed".
55 subroutine domain_read_1d(fileobj, variable_name, vdata, unlim_dim_level, &
56  corner, edge_lengths)
57 
58  type(fmsnetcdfdomainfile_t), intent(in) :: fileobj !< File object.
59  character(len=*), intent(in) :: variable_name !< Variable name.
60  class(*), dimension(:), intent(inout) :: vdata !< Data that will
61  !! be read out
62  !! to the netcdf file.
63  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
64  !! dimension.
65  integer, dimension(1), intent(in), optional :: corner !< Array of starting
66  !! indices describing
67  !! where the data
68  !! will be read to.
69  integer, dimension(1), intent(in), optional :: edge_lengths !< The number of
70  !! elements that
71  !! will be read
72  !! in each dimension.
73 
74  call netcdf_read_data(fileobj, variable_name, vdata, &
75  unlim_dim_level=unlim_dim_level, corner=corner, &
76  edge_lengths=edge_lengths, broadcast=.true.)
77 
78 end subroutine domain_read_1d
79 
80 
81 !> @brief I/O domain root reads in a domain decomposed variable at a
82 !! specific unlimited dimension level and scatters the data to the
83 !! rest of the ranks using its I/O compute domain indices. This
84 !! routine may only be used with variables that are "domain
85 !! decomposed".
86 subroutine domain_read_2d(fileobj, variable_name, vdata, unlim_dim_level, &
87  corner, edge_lengths)
88  type(fmsnetcdfdomainfile_t), intent(in) :: fileobj !< File object.
89  character(len=*), intent(in) :: variable_name !< Variable name.
90  class(*), contiguous, target, intent(inout) :: vdata(:,:) !< Data that will
91  !! be read out
92  !! to the netcdf file.
93  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
94  !! dimension.
95  integer, dimension(2), intent(in), optional :: corner !< Array of starting
96  !! indices describing
97  !! where the data
98  !! will be read to.
99  integer, dimension(2), intent(in), optional :: edge_lengths !< The number of
100  !! elements that
101  !! will be read
102  !! in each dimension.
103 
104  integer :: xdim_index !< The index of the variable that is the x dimension
105  integer :: ydim_index !< The index of the variable that is the y dimension
106  integer :: xpos !< The position of the x axis
107  integer :: ypos !< The position of the y axis
108  integer :: i !< For do loops
109  integer :: isd !< The starting x position of the data io_domain
110  integer :: isc !< The starting x position of the compute io_domain
111  integer :: xc_size !< The size of the x compute io_domain
112  integer :: yc_size !< The size of the y compute io_domain
113  integer :: jsd !< The ending x position of the data io_domain
114  integer :: jsc !< The ending y position of the compute io_domain
115  integer :: c(2) !< The corners of the data
116  integer :: e(2) !< The number of points (edges)
117  logical :: buffer_includes_halos !< .True. if vdata includes halo points
118  integer :: xgbegin !< Starting x index of global io domain
119  integer :: xgsize !< Size of global x io domain
120  integer :: ygbegin !< Starting y index of global io domain
121  integer :: ygsize !< Size of global y io domain
122  type(domain2d), pointer :: io_domain !< pointer to the io_domain
123 
124  !< The global data is only allocated by the io root PEs
125  integer(kind=i4_kind), dimension(:,:), allocatable :: buf_i4_kind_pe !< PES section of the data
126  integer(kind=i8_kind), dimension(:,:), allocatable :: buf_i8_kind_pe !< PES section of the data
127  real(kind=r4_kind), dimension(:,:), allocatable :: buf_r4_kind_pe !< PES section of the data
128  real(kind=r8_kind), dimension(:,:), allocatable :: buf_r8_kind_pe !< PES section of the data
129  integer(kind=i4_kind), dimension(:,:), allocatable :: buf_i4_kind !< Global section of the data
130  integer(kind=i8_kind), dimension(:,:), allocatable :: buf_i8_kind !< Global section of the data
131  real(kind=r4_kind), dimension(:,:), allocatable :: buf_r4_kind !< Global section of the data
132  real(kind=r8_kind), dimension(:,:), allocatable :: buf_r8_kind !< Global section of the data
133  class(*), dimension(:,:,:,:), pointer :: vdata_dummy !< Vdata remapped as 4D
134 
135  if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
136  xdim_index, ydim_index, xpos, ypos)) then
137  call netcdf_read_data(fileobj, variable_name, vdata, &
138  unlim_dim_level=unlim_dim_level, corner=corner, &
139  edge_lengths=edge_lengths, broadcast=.true.)
140  return
141  endif
142 
143  if (xdim_index .ne. 1 .or. ydim_index .ne. 2) then
144  ! This is a KLUDGE
145  ! mpp_scatter assumes that the variable is (x,y), if that is not the case it remaps the data
146  ! to a 4D array and calls domain_read_4d which does not use mpp_scatter yet
147  vdata_dummy(1:size(vdata,1),1:size(vdata,2), 1:1, 1:1) => vdata(:,:)
148  call domain_read_4d(fileobj, variable_name, vdata_dummy, unlim_dim_level)
149  return
150  endif
151  io_domain => mpp_get_io_domain(fileobj%domain)
152  call domain_offsets(size(vdata, xdim_index), size(vdata, ydim_index), fileobj%domain, &
153  xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
154  msg="file:"//trim(fileobj%path)//" and variable:"//trim(variable_name))
155  c(:) = 1
156  e(:) = shape(vdata)
157 
158  call mpp_get_global_domain(io_domain, xbegin=xgbegin, xsize=xgsize, position=xpos)
159  call mpp_get_global_domain(io_domain, ybegin=ygbegin, ysize=ygsize, position=ypos)
160 
161  !I/O root reads in the data and scatters it.
162  if (fileobj%is_root) then
163 
164  if (fileobj%adjust_indices) then
165  !< If the file is distributed, the file only contains the io global domain
166  c(xdim_index) = 1
167  c(ydim_index) = 1
168  else
169  !< If the file is not distributed read, the file contains the global domain,
170  !! so you only need to read the global io domain
171  c(xdim_index) = xgbegin
172  c(ydim_index) = ygbegin
173  endif
174 
175  e(xdim_index) = xgsize
176  e(ydim_index) = ygsize
177 
178  !Read in the global io domain
179  select type(vdata)
180  type is (integer(kind=i4_kind))
181  call allocate_array(buf_i4_kind, e)
182  call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
183  unlim_dim_level=unlim_dim_level, &
184  corner=c, edge_lengths=e, broadcast=.false.)
185  type is (integer(kind=i8_kind))
186  call allocate_array(buf_i8_kind, e)
187  call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
188  unlim_dim_level=unlim_dim_level, &
189  corner=c, edge_lengths=e, broadcast=.false.)
190  type is (real(kind=r4_kind))
191  call allocate_array(buf_r4_kind, e)
192  call netcdf_read_data(fileobj, variable_name, buf_r4_kind, &
193  unlim_dim_level=unlim_dim_level, &
194  corner=c, edge_lengths=e, broadcast=.false.)
195  type is (real(kind=r8_kind))
196  call allocate_array(buf_r8_kind, e)
197  call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
198  unlim_dim_level=unlim_dim_level, &
199  corner=c, edge_lengths=e, broadcast=.false.)
200  class default
201  call error("unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//" variable:"// &
202  & trim(variable_name))
203  end select
204 
205  endif
206 
207  c = 1
208  e = shape(vdata)
209 
210  if (buffer_includes_halos) then
211  !Adjust if the input buffer has room for halos.
212  c(xdim_index) = isc - isd + 1
213  c(ydim_index) = jsc - jsd + 1
214  else
215  c(xdim_index) = 1
216  c(ydim_index) = 1
217  endif
218 
219  e(xdim_index) = xc_size
220  e(ydim_index) = yc_size
221 
222  select type(vdata)
223  type is (integer(kind=i4_kind))
224  call allocate_array(buf_i4_kind_pe, e)
225  call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, fileobj%pelist, &
226  buf_i4_kind_pe, buf_i4_kind, fileobj%is_root)
227  call put_array_section(buf_i4_kind_pe, vdata, c, e)
228  deallocate(buf_i4_kind_pe)
229  type is (integer(kind=i8_kind))
230  call allocate_array(buf_i8_kind_pe, e)
231  call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, fileobj%pelist, &
232  buf_i8_kind_pe, buf_i8_kind, fileobj%is_root)
233  call put_array_section(buf_i8_kind_pe, vdata, c, e)
234  deallocate(buf_i8_kind_pe)
235  type is (real(kind=r4_kind))
236  call allocate_array(buf_r4_kind_pe, e)
237  call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, fileobj%pelist, &
238  buf_r4_kind_pe, buf_r4_kind, fileobj%is_root)
239  call put_array_section(buf_r4_kind_pe, vdata, c, e)
240  deallocate(buf_r4_kind_pe)
241  type is (real(kind=r8_kind))
242  call allocate_array(buf_r8_kind_pe, e)
243  call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, fileobj%pelist, &
244  buf_r8_kind_pe, buf_r8_kind, fileobj%is_root)
245  call put_array_section(buf_r8_kind_pe, vdata, c, e)
246  deallocate(buf_r8_kind_pe)
247  class default
248  call error("unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//" variable:"// &
249  & trim(variable_name))
250  end select
251 
252  if (fileobj%is_root) then
253  if (allocated(buf_i4_kind)) deallocate(buf_i4_kind)
254  if (allocated(buf_i8_kind)) deallocate(buf_i8_kind)
255  if (allocated(buf_r4_kind)) deallocate(buf_r4_kind)
256  if (allocated(buf_r8_kind)) deallocate(buf_r8_kind)
257  endif
258 end subroutine domain_read_2d
259 
260 
261 !> @brief I/O domain root reads in a domain decomposed variable at a
262 !! specific unlimited dimension level and scatters the data to the
263 !! rest of the ranks using its I/O compute domain indices. This
264 !! routine may only be used with variables that are "domain
265 !! decomposed".
266 subroutine domain_read_3d(fileobj, variable_name, vdata, unlim_dim_level, &
267  corner, edge_lengths)
268  type(fmsnetcdfdomainfile_t), intent(in) :: fileobj !< File object.
269  character(len=*), intent(in) :: variable_name !< Variable name.
270  class(*), contiguous, target, intent(inout) :: vdata(:,:,:) !< Data that will
271  !! be read out
272  !! to the netcdf file.
273  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
274  !! dimension.
275  integer, dimension(3), intent(in), optional :: corner !< Array of starting
276  !! indices describing
277  !! where the data
278  !! will be read to.
279  integer, dimension(3), intent(in), optional :: edge_lengths !< The number of
280  !! elements that
281  !! will be read
282  !! in each dimension.
283 
284  integer :: xdim_index !< The index of the variable that is the x dimension
285  integer :: ydim_index !< The index of the variable that is the y dimension
286  integer :: xpos !< The position of the x axis
287  integer :: ypos !< The position of the y axis
288  integer :: i !< For do loops
289  integer :: isd !< The starting x position of the data io_domain
290  integer :: isc !< The starting x position of the compute io_domain
291  integer :: xc_size !< The size of the x compute io_domain
292  integer :: yc_size !< The size of the y compute io_domain
293  integer :: jsd !< The ending x position of the data io_domain
294  integer :: jsc !< The ending y position of the compute io_domain
295  integer :: c(3) !< The corners of the data
296  integer :: e(3) !< The number of points (edges)
297  logical :: buffer_includes_halos !< .True. if vdata includes halo points
298  integer :: xgbegin !< Starting x index of global io domain
299  integer :: xgsize !< Size of global x io domain
300  integer :: ygbegin !< Starting y index of global io domain
301  integer :: ygsize !< Size of global y io domain
302  type(domain2d), pointer :: io_domain !< pointer to the io_domain
303 
304  !< The global data is only allocated by the io root PEs
305  integer(kind=i4_kind), dimension(:,:,:), allocatable :: buf_i4_kind_pe !< PES section of the data
306  integer(kind=i8_kind), dimension(:,:,:), allocatable :: buf_i8_kind_pe !< PES section of the data
307  real(kind=r4_kind), dimension(:,:,:), allocatable :: buf_r4_kind_pe !< PES section of the data
308  real(kind=r8_kind), dimension(:,:,:), allocatable :: buf_r8_kind_pe !< PES section of the data
309  integer(kind=i4_kind), dimension(:,:,:), allocatable :: buf_i4_kind !< Global section of the data
310  integer(kind=i8_kind), dimension(:,:,:), allocatable :: buf_i8_kind !< Global section of the data
311  real(kind=r4_kind), dimension(:,:,:), allocatable :: buf_r4_kind !< Global section of the data
312  real(kind=r8_kind), dimension(:,:,:), allocatable :: buf_r8_kind !< Global section of the data
313  class(*), dimension(:,:,:,:), pointer :: vdata_dummy !< Vdata remapped as 4D
314 
315  if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
316  xdim_index, ydim_index, xpos, ypos)) then
317  call netcdf_read_data(fileobj, variable_name, vdata, &
318  unlim_dim_level=unlim_dim_level, corner=corner, &
319  edge_lengths=edge_lengths, broadcast=.true.)
320  return
321  endif
322 
323  if (xdim_index .ne. 1 .or. ydim_index .ne. 2) then
324  ! This is a KLUDGE
325  ! mpp_scatter assumes that the variable is (x,y), if that is not the case it remaps the data
326  ! to a 4D array and calls domain_read_4d which does not use mpp_scatter yet
327  vdata_dummy(1:size(vdata,1),1:size(vdata,2), 1:size(vdata,3), 1:1) => vdata(:,:,:)
328  call domain_read_4d(fileobj, variable_name, vdata_dummy, unlim_dim_level)
329  return
330  endif
331  io_domain => mpp_get_io_domain(fileobj%domain)
332  call domain_offsets(size(vdata, xdim_index), size(vdata, ydim_index), fileobj%domain, &
333  xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
334  msg="file:"//trim(fileobj%path)//" and variable:"//trim(variable_name))
335  c(:) = 1
336  if (present(corner)) c = corner
337 
338  e(:) = shape(vdata)
339  if (present(edge_lengths)) e = edge_lengths
340 
341  call mpp_get_global_domain(io_domain, xbegin=xgbegin, xsize=xgsize, position=xpos)
342  call mpp_get_global_domain(io_domain, ybegin=ygbegin, ysize=ygsize, position=ypos)
343 
344  !I/O root reads in the data and scatters it.
345  if (fileobj%is_root) then
346 
347  if (fileobj%adjust_indices) then
348  !< If the file is distributed, the file only contains the io global domain
349  c(xdim_index) = 1
350  c(ydim_index) = 1
351  else
352  !< If the file is not distributed read, the file contains the global domain,
353  !! so you only need to read the global io domain
354  c(xdim_index) = xgbegin
355  c(ydim_index) = ygbegin
356  endif
357 
358  e(xdim_index) = xgsize
359  e(ydim_index) = ygsize
360 
361  !Read in the global io domain
362  select type(vdata)
363  type is (integer(kind=i4_kind))
364  call allocate_array(buf_i4_kind, e)
365  call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
366  unlim_dim_level=unlim_dim_level, &
367  corner=c, edge_lengths=e, broadcast=.false.)
368  type is (integer(kind=i8_kind))
369  call allocate_array(buf_i8_kind, e)
370  call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
371  unlim_dim_level=unlim_dim_level, &
372  corner=c, edge_lengths=e, broadcast=.false.)
373  type is (real(kind=r4_kind))
374  call allocate_array(buf_r4_kind, e)
375  call netcdf_read_data(fileobj, variable_name, buf_r4_kind, &
376  unlim_dim_level=unlim_dim_level, &
377  corner=c, edge_lengths=e, broadcast=.false.)
378  type is (real(kind=r8_kind))
379  call allocate_array(buf_r8_kind, e)
380  call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
381  unlim_dim_level=unlim_dim_level, &
382  corner=c, edge_lengths=e, broadcast=.false.)
383  class default
384  call error("unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//" variable:"// &
385  & trim(variable_name))
386  end select
387 
388  endif
389 
390  c = 1
391  e = shape(vdata)
392 
393  if (buffer_includes_halos) then
394  !Adjust if the input buffer has room for halos.
395  c(xdim_index) = isc - isd + 1
396  c(ydim_index) = jsc - jsd + 1
397  else
398  c(xdim_index) = 1
399  c(ydim_index) = 1
400  endif
401 
402  e(xdim_index) = xc_size
403  e(ydim_index) = yc_size
404 
405  select type(vdata)
406  type is (integer(kind=i4_kind))
407  call allocate_array(buf_i4_kind_pe, e)
408  call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
409  buf_i4_kind_pe, buf_i4_kind, fileobj%is_root)
410  call put_array_section(buf_i4_kind_pe, vdata, c, e)
411  deallocate(buf_i4_kind_pe)
412  type is (integer(kind=i8_kind))
413  call allocate_array(buf_i8_kind_pe, e)
414  call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
415  buf_i8_kind_pe, buf_i8_kind, fileobj%is_root)
416  call put_array_section(buf_i8_kind_pe, vdata, c, e)
417  deallocate(buf_i8_kind_pe)
418  type is (real(kind=r4_kind))
419  call allocate_array(buf_r4_kind_pe, e)
420  call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
421  buf_r4_kind_pe, buf_r4_kind, fileobj%is_root)
422  call put_array_section(buf_r4_kind_pe, vdata, c, e)
423  deallocate(buf_r4_kind_pe)
424  type is (real(kind=r8_kind))
425  call allocate_array(buf_r8_kind_pe, e)
426  call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
427  buf_r8_kind_pe, buf_r8_kind, fileobj%is_root)
428  call put_array_section(buf_r8_kind_pe, vdata, c, e)
429  deallocate(buf_r8_kind_pe)
430  class default
431  call error("unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//" variable:"// &
432  & trim(variable_name))
433  end select
434 
435  if (fileobj%is_root) then
436  if (allocated(buf_i4_kind)) deallocate(buf_i4_kind)
437  if (allocated(buf_i8_kind)) deallocate(buf_i8_kind)
438  if (allocated(buf_r4_kind)) deallocate(buf_r4_kind)
439  if (allocated(buf_r8_kind)) deallocate(buf_r8_kind)
440  endif
441 
442 end subroutine domain_read_3d
443 
444 
445 !> @brief I/O domain root reads in a domain decomposed variable at a
446 !! specific unlimited dimension level and scatters the data to the
447 !! rest of the ranks using its I/O compute domain indices. This
448 !! routine may only be used with variables that are "domain
449 !! decomposed".
450 subroutine domain_read_4d(fileobj, variable_name, vdata, unlim_dim_level, &
451  corner, edge_lengths)
452 
453  type(fmsnetcdfdomainfile_t), intent(in) :: fileobj !< File object.
454  character(len=*), intent(in) :: variable_name !< Variable name.
455  class(*), dimension(:,:,:,:), intent(inout) :: vdata !< Data that will
456  !! be read out
457  !! to the netcdf file.
458  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
459  !! dimension.
460  integer, dimension(4), intent(in), optional :: corner !< Array of starting
461  !! indices describing
462  !! where the data
463  !! will be read to.
464  integer, dimension(4), intent(in), optional :: edge_lengths !< The number of
465  !! elements that
466  !! will be read
467  !! in each dimension.
468 
469  integer :: xdim_index
470  integer :: ydim_index
471  type(domain2d), pointer :: io_domain
472  integer :: xpos
473  integer :: ypos
474  integer :: i
475  integer :: isd
476  integer :: isc
477  integer :: xc_size
478  integer :: jsd
479  integer :: jsc
480  integer :: yc_size
481  integer, dimension(:), allocatable :: pe_isc
482  integer, dimension(:), allocatable :: pe_icsize
483  integer, dimension(:), allocatable :: pe_jsc
484  integer, dimension(:), allocatable :: pe_jcsize
485  integer, dimension(4) :: c
486  integer, dimension(4) :: e
487  integer(kind=i4_kind), dimension(:,:,:,:), allocatable :: buf_i4_kind
488  integer(kind=i8_kind), dimension(:,:,:,:), allocatable :: buf_i8_kind
489  real(kind=r4_kind), dimension(:,:,:,:), allocatable :: buf_r4_kind
490  real(kind=r8_kind), dimension(:,:,:,:), allocatable :: buf_r8_kind
491  logical :: buffer_includes_halos
492  integer :: xgmin !< Starting x index of global io domain
493  integer :: ygmin !< Starting y index of global io domain
494 
495  if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
496  xdim_index, ydim_index, xpos, ypos)) then
497  call netcdf_read_data(fileobj, variable_name, vdata, &
498  unlim_dim_level=unlim_dim_level, corner=corner, &
499  edge_lengths=edge_lengths, broadcast=.true.)
500  return
501  endif
502  io_domain => mpp_get_io_domain(fileobj%domain)
503  call domain_offsets(size(vdata, xdim_index), size(vdata, ydim_index), fileobj%domain, &
504  xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
505  msg="file:"//trim(fileobj%path)//" and variable:"//trim(variable_name))
506  c(:) = 1
507  e(:) = shape(vdata)
508  if (present(edge_lengths)) e = edge_lengths
509 
510  !I/O root reads in the data and scatters it.
511  if (fileobj%is_root) then
512  allocate(pe_isc(size(fileobj%pelist)))
513  allocate(pe_icsize(size(fileobj%pelist)))
514  allocate(pe_jsc(size(fileobj%pelist)))
515  allocate(pe_jcsize(size(fileobj%pelist)))
516  call mpp_get_compute_domains(io_domain, xbegin=pe_isc, xsize=pe_icsize, position=xpos)
517  call mpp_get_compute_domains(io_domain, ybegin=pe_jsc, ysize=pe_jcsize, position=ypos)
518  call mpp_get_global_domain(io_domain, xbegin=xgmin, position=xpos)
519  call mpp_get_global_domain(io_domain, ybegin=ygmin, position=ypos)
520  do i = 1, size(fileobj%pelist)
521  if (present(corner)) c = corner
522  c(xdim_index) = pe_isc(i)
523  c(ydim_index) = pe_jsc(i)
524  if (fileobj%adjust_indices) then
525  c(xdim_index) = c(xdim_index) - xgmin + 1
526  c(ydim_index) = c(ydim_index) - ygmin + 1
527  endif
528  e(xdim_index) = pe_icsize(i)
529  e(ydim_index) = pe_jcsize(i)
530  select type(vdata)
531  type is (integer(kind=i4_kind))
532  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
533  call allocate_array(buf_i4_kind, e)
534  call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
535  unlim_dim_level=unlim_dim_level, &
536  corner=c, edge_lengths=e, broadcast=.false.)
537  if (i .eq. 1) then
538  !Root rank stores data directly.
539  c = 1
540  if (buffer_includes_halos) then
541  !Adjust if the input buffer has room for halos.
542  c(xdim_index) = isc - isd + 1
543  c(ydim_index) = jsc - jsd + 1
544  endif
545  call put_array_section(buf_i4_kind, vdata, c, e)
546  else
547  !Send data to non-root ranks.
548  call mpp_send(buf_i4_kind, size(buf_i4_kind), fileobj%pelist(i))
549  call mpp_sync_self(check=event_send)
550  endif
551  deallocate(buf_i4_kind)
552  type is (integer(kind=i8_kind))
553  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
554  call allocate_array(buf_i8_kind, e)
555  call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
556  unlim_dim_level=unlim_dim_level, &
557  corner=c, edge_lengths=e, broadcast=.false.)
558  if (i .eq. 1) then
559  !Root rank stores data directly.
560  c = 1
561  if (buffer_includes_halos) then
562  !Adjust if the input buffer has room for halos.
563  c(xdim_index) = isc - isd + 1
564  c(ydim_index) = jsc - jsd + 1
565  endif
566  call put_array_section(buf_i8_kind, vdata, c, e)
567  else
568  !Send data to non-root ranks.
569  call mpp_send(buf_i8_kind, size(buf_i8_kind), fileobj%pelist(i))
570  call mpp_sync_self(check=event_send)
571  endif
572  deallocate(buf_i8_kind)
573  type is (real(kind=r4_kind))
574  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
575  call allocate_array(buf_r4_kind, e)
576  call netcdf_read_data(fileobj, variable_name, buf_r4_kind, &
577  unlim_dim_level=unlim_dim_level, &
578  corner=c, edge_lengths=e, broadcast=.false.)
579  if (i .eq. 1) then
580  !Root rank stores data directly.
581  c = 1
582  if (buffer_includes_halos) then
583  !Adjust if the input buffer has room for halos.
584  c(xdim_index) = isc - isd + 1
585  c(ydim_index) = jsc - jsd + 1
586  endif
587  call put_array_section(buf_r4_kind, vdata, c, e)
588  else
589  !Send data to non-root ranks.
590  call mpp_send(buf_r4_kind, size(buf_r4_kind), fileobj%pelist(i))
591  call mpp_sync_self(check=event_send)
592  endif
593  deallocate(buf_r4_kind)
594  type is (real(kind=r8_kind))
595  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
596  call allocate_array(buf_r8_kind, e)
597  call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
598  unlim_dim_level=unlim_dim_level, &
599  corner=c, edge_lengths=e, broadcast=.false.)
600  if (i .eq. 1) then
601  !Root rank stores data directly.
602  c = 1
603  if (buffer_includes_halos) then
604  !Adjust if the input buffer has room for halos.
605  c(xdim_index) = isc - isd + 1
606  c(ydim_index) = jsc - jsd + 1
607  endif
608  call put_array_section(buf_r8_kind, vdata, c, e)
609  else
610  !Send data to non-root ranks.
611  call mpp_send(buf_r8_kind, size(buf_r8_kind), fileobj%pelist(i))
612  call mpp_sync_self(check=event_send)
613  endif
614  deallocate(buf_r8_kind)
615  class default
616  call error("unsupported variable type: domain_read_4d: file: "//trim(fileobj%path)//" variable:"// &
617  & trim(variable_name))
618  end select
619  enddo
620  deallocate(pe_isc)
621  deallocate(pe_icsize)
622  deallocate(pe_jsc)
623  deallocate(pe_jcsize)
624  else
625  c = 1
626  if (buffer_includes_halos) then
627  c(xdim_index) = isc - isd + 1
628  c(ydim_index) = jsc - jsd + 1
629  endif
630  e(xdim_index) = xc_size
631  e(ydim_index) = yc_size
632  select type(vdata)
633  type is (integer(kind=i4_kind))
634  call allocate_array(buf_i4_kind, e)
635  call mpp_recv(buf_i4_kind, size(buf_i4_kind), fileobj%io_root, block=.true.)
636  call put_array_section(buf_i4_kind, vdata, c, e)
637  deallocate(buf_i4_kind)
638  type is (integer(kind=i8_kind))
639  call allocate_array(buf_i8_kind, e)
640  call mpp_recv(buf_i8_kind, size(buf_i8_kind), fileobj%io_root, block=.true.)
641  call put_array_section(buf_i8_kind, vdata, c, e)
642  deallocate(buf_i8_kind)
643  type is (real(kind=r4_kind))
644  call allocate_array(buf_r4_kind, e)
645  call mpp_recv(buf_r4_kind, size(buf_r4_kind), fileobj%io_root, block=.true.)
646  call put_array_section(buf_r4_kind, vdata, c, e)
647  deallocate(buf_r4_kind)
648  type is (real(kind=r8_kind))
649  call allocate_array(buf_r8_kind, e)
650  call mpp_recv(buf_r8_kind, size(buf_r8_kind), fileobj%io_root, block=.true.)
651  call put_array_section(buf_r8_kind, vdata, c, e)
652  deallocate(buf_r8_kind)
653  class default
654  call error("unsupported variable type: domain_read_4d: file: "//trim(fileobj%path)//" variable:"// &
655  & trim(variable_name))
656  end select
657  endif
658 end subroutine domain_read_4d
659 
660 
661 !> @brief I/O domain root reads in a domain decomposed variable at a
662 !! specific unlimited dimension level and scatters the data to the
663 !! rest of the ranks using its I/O compute domain indices. This
664 !! routine may only be used with variables that are "domain
665 !! decomposed".
666 subroutine domain_read_5d(fileobj, variable_name, vdata, unlim_dim_level, &
667  corner, edge_lengths)
668 
669  type(fmsnetcdfdomainfile_t), intent(in) :: fileobj !< File object.
670  character(len=*), intent(in) :: variable_name !< Variable name.
671  class(*), dimension(:,:,:,:,:), intent(inout) :: vdata !< Data that will
672  !! be read out
673  !! to the netcdf file.
674  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
675  !! dimension.
676  integer, dimension(5), intent(in), optional :: corner !< Array of starting
677  !! indices describing
678  !! where the data
679  !! will be read to.
680  integer, dimension(5), intent(in), optional :: edge_lengths !< The number of
681  !! elements that
682  !! will be read
683  !! in each dimension.
684 
685  integer :: xdim_index
686  integer :: ydim_index
687  type(domain2d), pointer :: io_domain
688  integer :: xpos
689  integer :: ypos
690  integer :: i
691  integer :: isd
692  integer :: isc
693  integer :: xc_size
694  integer :: jsd
695  integer :: jsc
696  integer :: yc_size
697  integer, dimension(:), allocatable :: pe_isc
698  integer, dimension(:), allocatable :: pe_icsize
699  integer, dimension(:), allocatable :: pe_jsc
700  integer, dimension(:), allocatable :: pe_jcsize
701  integer, dimension(5) :: c
702  integer, dimension(5) :: e
703  integer(kind=i4_kind), dimension(:,:,:,:,:), allocatable :: buf_i4_kind
704  integer(kind=i8_kind), dimension(:,:,:,:,:), allocatable :: buf_i8_kind
705  real(kind=r4_kind), dimension(:,:,:,:,:), allocatable :: buf_r4_kind
706  real(kind=r8_kind), dimension(:,:,:,:,:), allocatable :: buf_r8_kind
707  logical :: buffer_includes_halos
708  integer :: xgmin !< Starting x index of global io domain
709  integer :: ygmin !< Starting y index of global io domain
710 
711  if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
712  xdim_index, ydim_index, xpos, ypos)) then
713  call netcdf_read_data(fileobj, variable_name, vdata, &
714  unlim_dim_level=unlim_dim_level, corner=corner, &
715  edge_lengths=edge_lengths, broadcast=.true.)
716  return
717  endif
718  io_domain => mpp_get_io_domain(fileobj%domain)
719  call domain_offsets(size(vdata, xdim_index), size(vdata, ydim_index), fileobj%domain, &
720  xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
721  msg="file:"//trim(fileobj%path)//" and variable:"//trim(variable_name))
722  c(:) = 1
723  e(:) = shape(vdata)
724  if (present(edge_lengths)) e = edge_lengths
725 
726  !I/O root reads in the data and scatters it.
727  if (fileobj%is_root) then
728  allocate(pe_isc(size(fileobj%pelist)))
729  allocate(pe_icsize(size(fileobj%pelist)))
730  allocate(pe_jsc(size(fileobj%pelist)))
731  allocate(pe_jcsize(size(fileobj%pelist)))
732  call mpp_get_compute_domains(io_domain, xbegin=pe_isc, xsize=pe_icsize, position=xpos)
733  call mpp_get_compute_domains(io_domain, ybegin=pe_jsc, ysize=pe_jcsize, position=ypos)
734  call mpp_get_global_domain(io_domain, xbegin=xgmin, position=xpos)
735  call mpp_get_global_domain(io_domain, ybegin=ygmin, position=ypos)
736  do i = 1, size(fileobj%pelist)
737  !Calculate the indices of the domain-decomposed chunk relative to its position in the file.
738  if (present(corner)) c = corner
739  c(xdim_index) = pe_isc(i)
740  c(ydim_index) = pe_jsc(i)
741  if (fileobj%adjust_indices) then
742  c(xdim_index) = c(xdim_index) - xgmin + 1
743  c(ydim_index) = c(ydim_index) - ygmin + 1
744  endif
745  e(xdim_index) = pe_icsize(i)
746  e(ydim_index) = pe_jcsize(i)
747  select type(vdata)
748  type is (integer(kind=i4_kind))
749  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
750  call allocate_array(buf_i4_kind, e)
751  call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
752  unlim_dim_level=unlim_dim_level, &
753  corner=c, edge_lengths=e, broadcast=.false.)
754  if (i .eq. 1) then
755  !Root rank stores data directly. Re-adjust the indicies relative
756  !to the input buffer vdata.
757  c = 1
758  if (buffer_includes_halos) then
759  !Adjust if the input buffer has room for halos.
760  c(xdim_index) = isc - isd + 1
761  c(ydim_index) = jsc - jsd + 1
762  endif
763  call put_array_section(buf_i4_kind, vdata, c, e)
764  else
765  !Send data to non-root ranks.
766  call mpp_send(buf_i4_kind, size(buf_i4_kind), fileobj%pelist(i))
767  call mpp_sync_self(check=event_send)
768  endif
769  deallocate(buf_i4_kind)
770  type is (integer(kind=i8_kind))
771  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
772  call allocate_array(buf_i8_kind, e)
773  call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
774  unlim_dim_level=unlim_dim_level, &
775  corner=c, edge_lengths=e, broadcast=.false.)
776  if (i .eq. 1) then
777  !Root rank stores data directly.
778  c = 1
779  if (buffer_includes_halos) then
780  !Adjust if the input buffer has room for halos.
781  c(xdim_index) = isc - isd + 1
782  c(ydim_index) = jsc - jsd + 1
783  endif
784  call put_array_section(buf_i8_kind, vdata, c, e)
785  else
786  !Send data to non-root ranks.
787  call mpp_send(buf_i8_kind, size(buf_i8_kind), fileobj%pelist(i))
788  call mpp_sync_self(check=event_send)
789  endif
790  deallocate(buf_i8_kind)
791  type is (real(kind=r4_kind))
792  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
793  call allocate_array(buf_r4_kind, e)
794  call netcdf_read_data(fileobj, variable_name, buf_r4_kind, &
795  unlim_dim_level=unlim_dim_level, &
796  corner=c, edge_lengths=e, broadcast=.false.)
797  if (i .eq. 1) then
798  !Root rank stores data directly.
799  c = 1
800  if (buffer_includes_halos) then
801  !Adjust if the input buffer has room for halos.
802  c(xdim_index) = isc - isd + 1
803  c(ydim_index) = jsc - jsd + 1
804  endif
805  call put_array_section(buf_r4_kind, vdata, c, e)
806  else
807  !Send data to non-root ranks.
808  call mpp_send(buf_r4_kind, size(buf_r4_kind), fileobj%pelist(i))
809  call mpp_sync_self(check=event_send)
810  endif
811  deallocate(buf_r4_kind)
812  type is (real(kind=r8_kind))
813  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
814  call allocate_array(buf_r8_kind, e)
815  call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
816  unlim_dim_level=unlim_dim_level, &
817  corner=c, edge_lengths=e, broadcast=.false.)
818  if (i .eq. 1) then
819  !Root rank stores data directly.
820  c = 1
821  if (buffer_includes_halos) then
822  !Adjust if the input buffer has room for halos.
823  c(xdim_index) = isc - isd + 1
824  c(ydim_index) = jsc - jsd + 1
825  endif
826  call put_array_section(buf_r8_kind, vdata, c, e)
827  else
828  !Send data to non-root ranks.
829  call mpp_send(buf_r8_kind, size(buf_r8_kind), fileobj%pelist(i))
830  call mpp_sync_self(check=event_send)
831  endif
832  deallocate(buf_r8_kind)
833  class default
834  call error("unsupported variable type: domain_read_5d: file: "//trim(fileobj%path)//" variable:"// &
835  & trim(variable_name))
836  end select
837  enddo
838  deallocate(pe_isc)
839  deallocate(pe_icsize)
840  deallocate(pe_jsc)
841  deallocate(pe_jcsize)
842  else
843  c = 1
844  if (buffer_includes_halos) then
845  c(xdim_index) = isc - isd + 1
846  c(ydim_index) = jsc - jsd + 1
847  endif
848  e(xdim_index) = xc_size
849  e(ydim_index) = yc_size
850  select type(vdata)
851  type is (integer(kind=i4_kind))
852  call allocate_array(buf_i4_kind, e)
853  call mpp_recv(buf_i4_kind, size(buf_i4_kind), fileobj%io_root, block=.true.)
854  call put_array_section(buf_i4_kind, vdata, c, e)
855  deallocate(buf_i4_kind)
856  type is (integer(kind=i8_kind))
857  call allocate_array(buf_i8_kind, e)
858  call mpp_recv(buf_i8_kind, size(buf_i8_kind), fileobj%io_root, block=.true.)
859  call put_array_section(buf_i8_kind, vdata, c, e)
860  deallocate(buf_i8_kind)
861  type is (real(kind=r4_kind))
862  call allocate_array(buf_r4_kind, e)
863  call mpp_recv(buf_r4_kind, size(buf_r4_kind), fileobj%io_root, block=.true.)
864  call put_array_section(buf_r4_kind, vdata, c, e)
865  deallocate(buf_r4_kind)
866  type is (real(kind=r8_kind))
867  call allocate_array(buf_r8_kind, e)
868  call mpp_recv(buf_r8_kind, size(buf_r8_kind), fileobj%io_root, block=.true.)
869  call put_array_section(buf_r8_kind, vdata, c, e)
870  deallocate(buf_r8_kind)
871  class default
872  call error("unsupported variable type: domain_read_5d: file: "//trim(fileobj%path)//" variable:"// &
873  & trim(variable_name))
874  end select
875  endif
876 end subroutine domain_read_5d
877 !> @}
subroutine domain_read_0d(fileobj, variable_name, vdata, unlim_dim_level, corner)
I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and sca...
Definition: domain_read.inc:30
subroutine domain_read_3d(fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and sca...
subroutine domain_read_1d(fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and sca...
Definition: domain_read.inc:57
subroutine domain_read_2d(fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and sca...
Definition: domain_read.inc:88
subroutine domain_read_4d(fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and sca...
subroutine domain_read_5d(fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and sca...
type(domain2d) function, pointer mpp_get_io_domain(domain)
Set user stack size.
subroutine mpp_sync_self(pelist, check, request, msg_size, msg_type)
This is to check if current PE's outstanding puts are complete but we can't use shmem_fence because w...