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