FMS  2024.03
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  e(:) = shape(vdata)
338 
339  call mpp_get_global_domain(io_domain, xbegin=xgbegin, xsize=xgsize, position=xpos)
340  call mpp_get_global_domain(io_domain, ybegin=ygbegin, ysize=ygsize, position=ypos)
341 
342  !I/O root reads in the data and scatters it.
343  if (fileobj%is_root) then
344 
345  if (fileobj%adjust_indices) then
346  !< If the file is distributed, the file only contains the io global domain
347  c(xdim_index) = 1
348  c(ydim_index) = 1
349  else
350  !< If the file is not distributed read, the file contains the global domain,
351  !! so you only need to read the global io domain
352  c(xdim_index) = xgbegin
353  c(ydim_index) = ygbegin
354  endif
355 
356  e(xdim_index) = xgsize
357  e(ydim_index) = ygsize
358 
359  !Read in the global io domain
360  select type(vdata)
361  type is (integer(kind=i4_kind))
362  call allocate_array(buf_i4_kind, e)
363  call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
364  unlim_dim_level=unlim_dim_level, &
365  corner=c, edge_lengths=e, broadcast=.false.)
366  type is (integer(kind=i8_kind))
367  call allocate_array(buf_i8_kind, e)
368  call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
369  unlim_dim_level=unlim_dim_level, &
370  corner=c, edge_lengths=e, broadcast=.false.)
371  type is (real(kind=r4_kind))
372  call allocate_array(buf_r4_kind, e)
373  call netcdf_read_data(fileobj, variable_name, buf_r4_kind, &
374  unlim_dim_level=unlim_dim_level, &
375  corner=c, edge_lengths=e, broadcast=.false.)
376  type is (real(kind=r8_kind))
377  call allocate_array(buf_r8_kind, e)
378  call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
379  unlim_dim_level=unlim_dim_level, &
380  corner=c, edge_lengths=e, broadcast=.false.)
381  class default
382  call error("unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//" variable:"// &
383  & trim(variable_name))
384  end select
385 
386  endif
387 
388  c = 1
389  e = shape(vdata)
390 
391  if (buffer_includes_halos) then
392  !Adjust if the input buffer has room for halos.
393  c(xdim_index) = isc - isd + 1
394  c(ydim_index) = jsc - jsd + 1
395  else
396  c(xdim_index) = 1
397  c(ydim_index) = 1
398  endif
399 
400  e(xdim_index) = xc_size
401  e(ydim_index) = yc_size
402 
403  select type(vdata)
404  type is (integer(kind=i4_kind))
405  call allocate_array(buf_i4_kind_pe, e)
406  call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
407  buf_i4_kind_pe, buf_i4_kind, fileobj%is_root)
408  call put_array_section(buf_i4_kind_pe, vdata, c, e)
409  deallocate(buf_i4_kind_pe)
410  type is (integer(kind=i8_kind))
411  call allocate_array(buf_i8_kind_pe, e)
412  call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
413  buf_i8_kind_pe, buf_i8_kind, fileobj%is_root)
414  call put_array_section(buf_i8_kind_pe, vdata, c, e)
415  deallocate(buf_i8_kind_pe)
416  type is (real(kind=r4_kind))
417  call allocate_array(buf_r4_kind_pe, e)
418  call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
419  buf_r4_kind_pe, buf_r4_kind, fileobj%is_root)
420  call put_array_section(buf_r4_kind_pe, vdata, c, e)
421  deallocate(buf_r4_kind_pe)
422  type is (real(kind=r8_kind))
423  call allocate_array(buf_r8_kind_pe, e)
424  call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
425  buf_r8_kind_pe, buf_r8_kind, fileobj%is_root)
426  call put_array_section(buf_r8_kind_pe, vdata, c, e)
427  deallocate(buf_r8_kind_pe)
428  class default
429  call error("unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//" variable:"// &
430  & trim(variable_name))
431  end select
432 
433  if (fileobj%is_root) then
434  if (allocated(buf_i4_kind)) deallocate(buf_i4_kind)
435  if (allocated(buf_i8_kind)) deallocate(buf_i8_kind)
436  if (allocated(buf_r4_kind)) deallocate(buf_r4_kind)
437  if (allocated(buf_r8_kind)) deallocate(buf_r8_kind)
438  endif
439 
440 end subroutine domain_read_3d
441 
442 
443 !> @brief I/O domain root reads in a domain decomposed variable at a
444 !! specific unlimited dimension level and scatters the data to the
445 !! rest of the ranks using its I/O compute domain indices. This
446 !! routine may only be used with variables that are "domain
447 !! decomposed".
448 subroutine domain_read_4d(fileobj, variable_name, vdata, unlim_dim_level, &
449  corner, edge_lengths)
450 
451  type(fmsnetcdfdomainfile_t), intent(in) :: fileobj !< File object.
452  character(len=*), intent(in) :: variable_name !< Variable name.
453  class(*), dimension(:,:,:,:), intent(inout) :: vdata !< Data that will
454  !! be read out
455  !! to the netcdf file.
456  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
457  !! dimension.
458  integer, dimension(4), intent(in), optional :: corner !< Array of starting
459  !! indices describing
460  !! where the data
461  !! will be read to.
462  integer, dimension(4), intent(in), optional :: edge_lengths !< The number of
463  !! elements that
464  !! will be read
465  !! in each dimension.
466 
467  integer :: xdim_index
468  integer :: ydim_index
469  type(domain2d), pointer :: io_domain
470  integer :: xpos
471  integer :: ypos
472  integer :: i
473  integer :: isd
474  integer :: isc
475  integer :: xc_size
476  integer :: jsd
477  integer :: jsc
478  integer :: yc_size
479  integer, dimension(:), allocatable :: pe_isc
480  integer, dimension(:), allocatable :: pe_icsize
481  integer, dimension(:), allocatable :: pe_jsc
482  integer, dimension(:), allocatable :: pe_jcsize
483  integer, dimension(4) :: c
484  integer, dimension(4) :: e
485  integer(kind=i4_kind), dimension(:,:,:,:), allocatable :: buf_i4_kind
486  integer(kind=i8_kind), dimension(:,:,:,:), allocatable :: buf_i8_kind
487  real(kind=r4_kind), dimension(:,:,:,:), allocatable :: buf_r4_kind
488  real(kind=r8_kind), dimension(:,:,:,:), allocatable :: buf_r8_kind
489  logical :: buffer_includes_halos
490  integer :: xgmin !< Starting x index of global io domain
491  integer :: ygmin !< Starting y index of global io domain
492 
493  if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
494  xdim_index, ydim_index, xpos, ypos)) then
495  call netcdf_read_data(fileobj, variable_name, vdata, &
496  unlim_dim_level=unlim_dim_level, corner=corner, &
497  edge_lengths=edge_lengths, broadcast=.true.)
498  return
499  endif
500  io_domain => mpp_get_io_domain(fileobj%domain)
501  call domain_offsets(size(vdata, xdim_index), size(vdata, ydim_index), fileobj%domain, &
502  xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
503  msg="file:"//trim(fileobj%path)//" and variable:"//trim(variable_name))
504  c(:) = 1
505  e(:) = shape(vdata)
506 
507  !I/O root reads in the data and scatters it.
508  if (fileobj%is_root) then
509  allocate(pe_isc(size(fileobj%pelist)))
510  allocate(pe_icsize(size(fileobj%pelist)))
511  allocate(pe_jsc(size(fileobj%pelist)))
512  allocate(pe_jcsize(size(fileobj%pelist)))
513  call mpp_get_compute_domains(io_domain, xbegin=pe_isc, xsize=pe_icsize, position=xpos)
514  call mpp_get_compute_domains(io_domain, ybegin=pe_jsc, ysize=pe_jcsize, position=ypos)
515  call mpp_get_global_domain(io_domain, xbegin=xgmin, position=xpos)
516  call mpp_get_global_domain(io_domain, ybegin=ygmin, position=ypos)
517  do i = 1, size(fileobj%pelist)
518  c(xdim_index) = pe_isc(i)
519  c(ydim_index) = pe_jsc(i)
520  if (fileobj%adjust_indices) then
521  c(xdim_index) = c(xdim_index) - xgmin + 1
522  c(ydim_index) = c(ydim_index) - ygmin + 1
523  endif
524  e(xdim_index) = pe_icsize(i)
525  e(ydim_index) = pe_jcsize(i)
526  select type(vdata)
527  type is (integer(kind=i4_kind))
528  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
529  call allocate_array(buf_i4_kind, e)
530  call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
531  unlim_dim_level=unlim_dim_level, &
532  corner=c, edge_lengths=e, broadcast=.false.)
533  if (i .eq. 1) then
534  !Root rank stores data directly.
535  if (buffer_includes_halos) then
536  !Adjust if the input buffer has room for halos.
537  c(xdim_index) = isc - isd + 1
538  c(ydim_index) = jsc - jsd + 1
539  else
540  c(xdim_index) = 1
541  c(ydim_index) = 1
542  endif
543  call put_array_section(buf_i4_kind, vdata, c, e)
544  else
545  !Send data to non-root ranks.
546  call mpp_send(buf_i4_kind, size(buf_i4_kind), fileobj%pelist(i))
547  call mpp_sync_self(check=event_send)
548  endif
549  deallocate(buf_i4_kind)
550  type is (integer(kind=i8_kind))
551  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
552  call allocate_array(buf_i8_kind, e)
553  call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
554  unlim_dim_level=unlim_dim_level, &
555  corner=c, edge_lengths=e, broadcast=.false.)
556  if (i .eq. 1) then
557  !Root rank stores data directly.
558  if (buffer_includes_halos) then
559  !Adjust if the input buffer has room for halos.
560  c(xdim_index) = isc - isd + 1
561  c(ydim_index) = jsc - jsd + 1
562  else
563  c(xdim_index) = 1
564  c(ydim_index) = 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  if (buffer_includes_halos) then
582  !Adjust if the input buffer has room for halos.
583  c(xdim_index) = isc - isd + 1
584  c(ydim_index) = jsc - jsd + 1
585  else
586  c(xdim_index) = 1
587  c(ydim_index) = 1
588  endif
589  call put_array_section(buf_r4_kind, vdata, c, e)
590  else
591  !Send data to non-root ranks.
592  call mpp_send(buf_r4_kind, size(buf_r4_kind), fileobj%pelist(i))
593  call mpp_sync_self(check=event_send)
594  endif
595  deallocate(buf_r4_kind)
596  type is (real(kind=r8_kind))
597  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
598  call allocate_array(buf_r8_kind, e)
599  call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
600  unlim_dim_level=unlim_dim_level, &
601  corner=c, edge_lengths=e, broadcast=.false.)
602  if (i .eq. 1) then
603  !Root rank stores data directly.
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  else
609  c(xdim_index) = 1
610  c(ydim_index) = 1
611  endif
612  call put_array_section(buf_r8_kind, vdata, c, e)
613  else
614  !Send data to non-root ranks.
615  call mpp_send(buf_r8_kind, size(buf_r8_kind), fileobj%pelist(i))
616  call mpp_sync_self(check=event_send)
617  endif
618  deallocate(buf_r8_kind)
619  class default
620  call error("unsupported variable type: domain_read_4d: file: "//trim(fileobj%path)//" variable:"// &
621  & trim(variable_name))
622  end select
623  enddo
624  deallocate(pe_isc)
625  deallocate(pe_icsize)
626  deallocate(pe_jsc)
627  deallocate(pe_jcsize)
628  else
629  if (buffer_includes_halos) then
630  c(xdim_index) = isc - isd + 1
631  c(ydim_index) = jsc - jsd + 1
632  endif
633  e(xdim_index) = xc_size
634  e(ydim_index) = yc_size
635  select type(vdata)
636  type is (integer(kind=i4_kind))
637  call allocate_array(buf_i4_kind, e)
638  call mpp_recv(buf_i4_kind, size(buf_i4_kind), fileobj%io_root, block=.true.)
639  call put_array_section(buf_i4_kind, vdata, c, e)
640  deallocate(buf_i4_kind)
641  type is (integer(kind=i8_kind))
642  call allocate_array(buf_i8_kind, e)
643  call mpp_recv(buf_i8_kind, size(buf_i8_kind), fileobj%io_root, block=.true.)
644  call put_array_section(buf_i8_kind, vdata, c, e)
645  deallocate(buf_i8_kind)
646  type is (real(kind=r4_kind))
647  call allocate_array(buf_r4_kind, e)
648  call mpp_recv(buf_r4_kind, size(buf_r4_kind), fileobj%io_root, block=.true.)
649  call put_array_section(buf_r4_kind, vdata, c, e)
650  deallocate(buf_r4_kind)
651  type is (real(kind=r8_kind))
652  call allocate_array(buf_r8_kind, e)
653  call mpp_recv(buf_r8_kind, size(buf_r8_kind), fileobj%io_root, block=.true.)
654  call put_array_section(buf_r8_kind, vdata, c, e)
655  deallocate(buf_r8_kind)
656  class default
657  call error("unsupported variable type: domain_read_4d: file: "//trim(fileobj%path)//" variable:"// &
658  & trim(variable_name))
659  end select
660  endif
661 end subroutine domain_read_4d
662 
663 
664 !> @brief I/O domain root reads in a domain decomposed variable at a
665 !! specific unlimited dimension level and scatters the data to the
666 !! rest of the ranks using its I/O compute domain indices. This
667 !! routine may only be used with variables that are "domain
668 !! decomposed".
669 subroutine domain_read_5d(fileobj, variable_name, vdata, unlim_dim_level, &
670  corner, edge_lengths)
671 
672  type(fmsnetcdfdomainfile_t), intent(in) :: fileobj !< File object.
673  character(len=*), intent(in) :: variable_name !< Variable name.
674  class(*), dimension(:,:,:,:,:), intent(inout) :: vdata !< Data that will
675  !! be read out
676  !! to the netcdf file.
677  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
678  !! dimension.
679  integer, dimension(5), intent(in), optional :: corner !< Array of starting
680  !! indices describing
681  !! where the data
682  !! will be read to.
683  integer, dimension(5), intent(in), optional :: edge_lengths !< The number of
684  !! elements that
685  !! will be read
686  !! in each dimension.
687 
688  integer :: xdim_index
689  integer :: ydim_index
690  type(domain2d), pointer :: io_domain
691  integer :: xpos
692  integer :: ypos
693  integer :: i
694  integer :: isd
695  integer :: isc
696  integer :: xc_size
697  integer :: jsd
698  integer :: jsc
699  integer :: yc_size
700  integer, dimension(:), allocatable :: pe_isc
701  integer, dimension(:), allocatable :: pe_icsize
702  integer, dimension(:), allocatable :: pe_jsc
703  integer, dimension(:), allocatable :: pe_jcsize
704  integer, dimension(5) :: c
705  integer, dimension(5) :: e
706  integer(kind=i4_kind), dimension(:,:,:,:,:), allocatable :: buf_i4_kind
707  integer(kind=i8_kind), dimension(:,:,:,:,:), allocatable :: buf_i8_kind
708  real(kind=r4_kind), dimension(:,:,:,:,:), allocatable :: buf_r4_kind
709  real(kind=r8_kind), dimension(:,:,:,:,:), allocatable :: buf_r8_kind
710  logical :: buffer_includes_halos
711  integer :: xgmin !< Starting x index of global io domain
712  integer :: ygmin !< Starting y index of global io domain
713 
714  if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
715  xdim_index, ydim_index, xpos, ypos)) then
716  call netcdf_read_data(fileobj, variable_name, vdata, &
717  unlim_dim_level=unlim_dim_level, corner=corner, &
718  edge_lengths=edge_lengths, broadcast=.true.)
719  return
720  endif
721  io_domain => mpp_get_io_domain(fileobj%domain)
722  call domain_offsets(size(vdata, xdim_index), size(vdata, ydim_index), fileobj%domain, &
723  xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
724  msg="file:"//trim(fileobj%path)//" and variable:"//trim(variable_name))
725  c(:) = 1
726  e(:) = shape(vdata)
727 
728  !I/O root reads in the data and scatters it.
729  if (fileobj%is_root) then
730  allocate(pe_isc(size(fileobj%pelist)))
731  allocate(pe_icsize(size(fileobj%pelist)))
732  allocate(pe_jsc(size(fileobj%pelist)))
733  allocate(pe_jcsize(size(fileobj%pelist)))
734  call mpp_get_compute_domains(io_domain, xbegin=pe_isc, xsize=pe_icsize, position=xpos)
735  call mpp_get_compute_domains(io_domain, ybegin=pe_jsc, ysize=pe_jcsize, position=ypos)
736  call mpp_get_global_domain(io_domain, xbegin=xgmin, position=xpos)
737  call mpp_get_global_domain(io_domain, ybegin=ygmin, position=ypos)
738  do i = 1, size(fileobj%pelist)
739  !Calculate the indices of the domain-decomposed chunk relative to its position in the file.
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  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  else
763  c(xdim_index) = 1
764  c(ydim_index) = 1
765  endif
766  call put_array_section(buf_i4_kind, vdata, c, e)
767  else
768  !Send data to non-root ranks.
769  call mpp_send(buf_i4_kind, size(buf_i4_kind), fileobj%pelist(i))
770  call mpp_sync_self(check=event_send)
771  endif
772  deallocate(buf_i4_kind)
773  type is (integer(kind=i8_kind))
774  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
775  call allocate_array(buf_i8_kind, e)
776  call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
777  unlim_dim_level=unlim_dim_level, &
778  corner=c, edge_lengths=e, broadcast=.false.)
779  if (i .eq. 1) then
780  !Root rank stores data directly.
781  if (buffer_includes_halos) then
782  !Adjust if the input buffer has room for halos.
783  c(xdim_index) = isc - isd + 1
784  c(ydim_index) = jsc - jsd + 1
785  else
786  c(xdim_index) = 1
787  c(ydim_index) = 1
788  endif
789  call put_array_section(buf_i8_kind, vdata, c, e)
790  else
791  !Send data to non-root ranks.
792  call mpp_send(buf_i8_kind, size(buf_i8_kind), fileobj%pelist(i))
793  call mpp_sync_self(check=event_send)
794  endif
795  deallocate(buf_i8_kind)
796  type is (real(kind=r4_kind))
797  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
798  call allocate_array(buf_r4_kind, e)
799  call netcdf_read_data(fileobj, variable_name, buf_r4_kind, &
800  unlim_dim_level=unlim_dim_level, &
801  corner=c, edge_lengths=e, broadcast=.false.)
802  if (i .eq. 1) then
803  !Root rank stores data directly.
804  if (buffer_includes_halos) then
805  !Adjust if the input buffer has room for halos.
806  c(xdim_index) = isc - isd + 1
807  c(ydim_index) = jsc - jsd + 1
808  else
809  c(xdim_index) = 1
810  c(ydim_index) = 1
811  endif
812  call put_array_section(buf_r4_kind, vdata, c, e)
813  else
814  !Send data to non-root ranks.
815  call mpp_send(buf_r4_kind, size(buf_r4_kind), fileobj%pelist(i))
816  call mpp_sync_self(check=event_send)
817  endif
818  deallocate(buf_r4_kind)
819  type is (real(kind=r8_kind))
820  !Read in the data for fileobj%pelist(i)'s portion of the compute domain.
821  call allocate_array(buf_r8_kind, e)
822  call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
823  unlim_dim_level=unlim_dim_level, &
824  corner=c, edge_lengths=e, broadcast=.false.)
825  if (i .eq. 1) then
826  !Root rank stores data directly.
827  if (buffer_includes_halos) then
828  !Adjust if the input buffer has room for halos.
829  c(xdim_index) = isc - isd + 1
830  c(ydim_index) = jsc - jsd + 1
831  else
832  c(xdim_index) = 1
833  c(ydim_index) = 1
834  endif
835  call put_array_section(buf_r8_kind, vdata, c, e)
836  else
837  !Send data to non-root ranks.
838  call mpp_send(buf_r8_kind, size(buf_r8_kind), fileobj%pelist(i))
839  call mpp_sync_self(check=event_send)
840  endif
841  deallocate(buf_r8_kind)
842  class default
843  call error("unsupported variable type: domain_read_5d: file: "//trim(fileobj%path)//" variable:"// &
844  & trim(variable_name))
845  end select
846  enddo
847  deallocate(pe_isc)
848  deallocate(pe_icsize)
849  deallocate(pe_jsc)
850  deallocate(pe_jcsize)
851  else
852  if (buffer_includes_halos) then
853  c(xdim_index) = isc - isd + 1
854  c(ydim_index) = jsc - jsd + 1
855  endif
856  e(xdim_index) = xc_size
857  e(ydim_index) = yc_size
858  select type(vdata)
859  type is (integer(kind=i4_kind))
860  call allocate_array(buf_i4_kind, e)
861  call mpp_recv(buf_i4_kind, size(buf_i4_kind), fileobj%io_root, block=.true.)
862  call put_array_section(buf_i4_kind, vdata, c, e)
863  deallocate(buf_i4_kind)
864  type is (integer(kind=i8_kind))
865  call allocate_array(buf_i8_kind, e)
866  call mpp_recv(buf_i8_kind, size(buf_i8_kind), fileobj%io_root, block=.true.)
867  call put_array_section(buf_i8_kind, vdata, c, e)
868  deallocate(buf_i8_kind)
869  type is (real(kind=r4_kind))
870  call allocate_array(buf_r4_kind, e)
871  call mpp_recv(buf_r4_kind, size(buf_r4_kind), fileobj%io_root, block=.true.)
872  call put_array_section(buf_r4_kind, vdata, c, e)
873  deallocate(buf_r4_kind)
874  type is (real(kind=r8_kind))
875  call allocate_array(buf_r8_kind, e)
876  call mpp_recv(buf_r8_kind, size(buf_r8_kind), fileobj%io_root, block=.true.)
877  call put_array_section(buf_r8_kind, vdata, c, e)
878  deallocate(buf_r8_kind)
879  class default
880  call error("unsupported variable type: domain_read_5d: file: "//trim(fileobj%path)//" variable:"// &
881  & trim(variable_name))
882  end select
883  endif
884 end subroutine domain_read_5d
885 !> @}
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...