FMS  2025.02
Flexible Modeling System
compressed_write.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 Compressed write routines for @ref write_data interface
21 
22 !> @addtogroup fms2_io_mod
23 !> @{
24 
25 !> @brief For variables without a compressed dimension, this routine simply wraps
26 !! netcdf_write data. If the variable does have a compressed axis, the I/O root
27 !! gathers the data from the rest of the ranks and then writes the combined data to
28 !! the netcdf file.
29 subroutine compressed_write_0d(fileobj, variable_name, cdata, unlim_dim_level, &
30  corner)
31 
32  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
33  character(len=*), intent(in) :: variable_name !< Variable name.
34  class(*), intent(in) :: cdata !< Compressed data that
35  !! will be gathered and
36  !! written to the
37  !! netcdf file.
38  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
39  !! dimension.
40  integer, intent(in), optional :: corner !< Array of starting
41  !! indices describing
42  !! where the data
43  !! will be written to.
44 
45  integer, dimension(2) :: compressed_dim_index
46 
47  compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
48  if (compressed_dim_index(1) .eq. dimension_not_found) then
49  call netcdf_write_data(fileobj, variable_name, cdata, &
50  unlim_dim_level=unlim_dim_level, corner=corner)
51  return
52  endif
53 end subroutine compressed_write_0d
54 
55 
56 !> @brief Wrapper to distinguish interfaces.
57 subroutine compressed_write_0d_wrap(fileobj, variable_name, cdata, unlim_dim_level, &
58  corner)
59 
60  type(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
61  character(len=*), intent(in) :: variable_name !< Variable name.
62  class(*), intent(in) :: cdata !< Compressed data that
63  !! will be gathered and
64  !! written to the
65  !! netcdf file.
66  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
67  !! dimension.
68  integer, intent(in), optional :: corner !< Array of starting
69  !! indices describing
70  !! where the data
71  !! will be written to.
72 
73  call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner)
74 end subroutine compressed_write_0d_wrap
75 
76 
77 !> @brief For variables without a compressed dimension, this routine simply wraps
78 !! netcdf_write data. If the variable does have a compressed axis, the I/O root
79 !! gathers the data from the rest of the ranks and then writes the combined data to
80 !! the netcdf file.
81 subroutine compressed_write_1d(fileobj, variable_name, cdata, unlim_dim_level, &
82  corner, edge_lengths)
83 
84  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
85  character(len=*), intent(in) :: variable_name !< Variable name.
86  class(*), dimension(:), intent(in) :: cdata !< Compressed data that
87  !! will be gathered and
88  !! written to the
89  !! netcdf file.
90  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
91  !! dimension.
92  integer, dimension(1), intent(in), optional :: corner !< Array of starting
93  !! indices describing
94  !! where the data
95  !! will be written to.
96  integer, dimension(1), intent(in), optional :: edge_lengths !< The number of
97  !! elements that
98  !! will be written
99  !! in each dimension.
100 
101  integer, dimension(2) :: compressed_dim_index !< index of the compressed dimension
102  !! compressed_dim_index(1) relative to cdata
103  !! compressed_dim_index(2) relative to the fileobj
104  integer, dimension(1) :: e !< "edges" number of points to read
105 
106  integer(kind=i4_kind), dimension(:), allocatable :: buf_i4_kind !< Global buffer of data
107  integer(kind=i8_kind), dimension(:), allocatable :: buf_i8_kind !< Global buffer of data
108  real (kind=r4_kind), dimension(:), allocatable :: buf_r4_kind !< Global buffer of data
109  real (kind=r8_kind), dimension(:), allocatable :: buf_r8_kind !< Global buffer of data
110 
111  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
112 
113  append_error_msg = "compressed_write_1d: file:"//trim(fileobj%path)//" variable:"//trim(variable_name)
114 
115  compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
116 
117  if (compressed_dim_index(1) .eq. dimension_not_found) then
118  call netcdf_write_data(fileobj, variable_name, cdata, &
119  unlim_dim_level=unlim_dim_level, corner=corner, &
120  edge_lengths=edge_lengths)
121  return
122  endif
123 
124  e(:) = shape(cdata)
125  !The root pe creates a buffer big enough to store the data:
126  if (fileobj%is_root) then
127  e(compressed_dim_index(1)) = sum(fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems)
128  select type(cdata)
129  type is (integer(kind=i4_kind))
130  call allocate_array(buf_i4_kind, e)
131  type is (integer(kind=i8_kind))
132  call allocate_array(buf_i8_kind, e)
133  type is (real(kind=r4_kind))
134  call allocate_array(buf_r4_kind, e)
135  type is (real(kind=r8_kind))
136  call allocate_array(buf_r8_kind, e)
137  class default
138  call error("unsupported variable type: "//trim(append_error_msg))
139  end select
140  else
141  select type(cdata)
142  type is (integer(kind=i4_kind))
143  allocate(buf_i4_kind(1))
144  type is (integer(kind=i8_kind))
145  allocate(buf_i8_kind(1))
146  type is (real(kind=r4_kind))
147  allocate(buf_r4_kind(1))
148  type is (real(kind=r8_kind))
149  allocate(buf_r8_kind(1))
150  class default
151  call error("unsupported variable type: "//trim(append_error_msg))
152  end select
153  endif
154 
155  !Gather the data onto the I/O root and write it out.
156  select type(cdata)
157  type is (integer(kind=i4_kind))
158  call mpp_gather(cdata, size(cdata), buf_i4_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
159  fileobj%pelist)
160  if (fileobj%is_root) call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
161  unlim_dim_level=unlim_dim_level)
162  type is (integer(kind=i8_kind))
163  call mpp_gather(cdata, size(cdata), buf_i8_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
164  fileobj%pelist)
165  if (fileobj%is_root) call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
166  unlim_dim_level=unlim_dim_level)
167  type is (real(kind=r4_kind))
168  call mpp_gather(cdata, size(cdata), buf_r4_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
169  fileobj%pelist)
170  if (fileobj%is_root) call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
171  unlim_dim_level=unlim_dim_level)
172  type is (real(kind=r8_kind))
173  call mpp_gather(cdata, size(cdata), buf_r8_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
174  fileobj%pelist)
175  if (fileobj%is_root) call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
176  unlim_dim_level=unlim_dim_level)
177  class default
178  call error("unsupported variable type: "//trim(append_error_msg))
179  end select
180  if (allocated(buf_i4_kind)) deallocate(buf_i4_kind)
181  if (allocated(buf_i8_kind)) deallocate(buf_i8_kind)
182  if (allocated(buf_r4_kind)) deallocate(buf_r4_kind)
183  if (allocated(buf_r8_kind)) deallocate(buf_r8_kind)
184 end subroutine compressed_write_1d
185 
186 
187 !> @brief For variables without a compressed dimension, this routine simply wraps
188 !! netcdf_write data. If the variable does have a compressed axis, the I/O root
189 !! gathers the data from the rest of the ranks and then writes the combined data to
190 !! the netcdf file.
191 subroutine compressed_write_2d(fileobj, variable_name, cdata, unlim_dim_level, &
192  corner, edge_lengths)
193 
194  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
195  character(len=*), intent(in) :: variable_name !< Variable name.
196  class(*), dimension(:,:), intent(in) :: cdata !< Compressed data that
197  !! will be gathered and
198  !! written to the
199  !! netcdf file.
200  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
201  !! dimension.
202  integer, dimension(2), intent(in), optional :: corner !< Array of starting
203  !! indices describing
204  !! where the data
205  !! will be written to.
206  integer, dimension(2), intent(in), optional :: edge_lengths !< The number of
207  !! elements that
208  !! will be written
209  !! in each dimension.
210 
211  integer, dimension(2) :: compressed_dim_index !< index of the compressed dimension
212  !! compressed_dim_index(1) relative to cdata
213  !! compressed_dim_index(2) relative to the fileobj
214  integer, dimension(2) :: c !! corners of the data to read
215  integer, dimension(2) :: e !< "edges" number of points to read
216 
217  integer(kind=i4_kind), dimension(:,:), allocatable :: buf_i4_kind !< Global buffer of data
218  integer(kind=i8_kind), dimension(:,:), allocatable :: buf_i8_kind !< Global buffer of data
219  real (kind=r4_kind), dimension(:,:), allocatable :: buf_r4_kind !< Global buffer of data
220  real (kind=r8_kind), dimension(:,:), allocatable :: buf_r8_kind !< Global buffer of data
221 
222  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
223 
224  integer :: index(1) !< index of the PE in the pelist
225 
226  integer :: is !< Starting index of the first dimension
227  integer :: ie !< Ending index of the first dimension
228  integer :: js !< Starting index of the second dimension
229  integer :: je !< Ending index of the second dimension
230 
231  append_error_msg = "compressed_write_2d: file:"//trim(fileobj%path)//" variable:"//trim(variable_name)
232 
233  compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
234  if (compressed_dim_index(1) .eq. dimension_not_found) then
235  call netcdf_write_data(fileobj, variable_name, cdata, &
236  unlim_dim_level=unlim_dim_level, corner=corner, &
237  edge_lengths=edge_lengths)
238  return
239  endif
240 
241  e(:) = shape(cdata)
242  !The root pe creates a buffer big enough to store the data:
243  if (fileobj%is_root) then
244  e(compressed_dim_index(1)) = sum(fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems)
245  select type(cdata)
246  type is (integer(kind=i4_kind))
247  call allocate_array(buf_i4_kind, e)
248  type is (integer(kind=i8_kind))
249  call allocate_array(buf_i8_kind, e)
250  type is (real(kind=r4_kind))
251  call allocate_array(buf_r4_kind, e)
252  type is (real(kind=r8_kind))
253  call allocate_array(buf_r8_kind, e)
254  class default
255  call error("unsupported variable type: "//trim(append_error_msg))
256  end select
257  endif
258 
259  c(:) = 1
260  index = findloc(fileobj%pelist, mpp_pe())
261 
262  c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(index(1))
263  e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(index(1))
264 
265  if (compressed_dim_index(1) .eq. 1) then
266  is = c(compressed_dim_index(1))
267  ie = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
268  js = 1
269  je = size(cdata,2)
270  else
271  js = c(compressed_dim_index(1))
272  je = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
273  is = 1
274  ie = size(cdata,2)
275  endif
276 
277  select type(cdata)
278  type is (integer(kind=i4_kind))
279  call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_i4_kind, fileobj%is_root)
280  if (fileobj%is_root) call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
281  unlim_dim_level=unlim_dim_level)
282  type is (integer(kind=i8_kind))
283  call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_i8_kind, fileobj%is_root)
284  if (fileobj%is_root) call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
285  unlim_dim_level=unlim_dim_level)
286  type is (real(kind=r4_kind))
287  call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_r4_kind, fileobj%is_root)
288  if (fileobj%is_root) call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
289  unlim_dim_level=unlim_dim_level)
290  type is (real(kind=r8_kind))
291  call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_r8_kind, fileobj%is_root)
292  if (fileobj%is_root) call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
293  unlim_dim_level=unlim_dim_level)
294  class default
295  call error("unsupported variable type: "//trim(append_error_msg))
296  end select
297 
298  if (fileobj%is_root) then
299  if (allocated(buf_i4_kind)) deallocate(buf_i4_kind)
300  if (allocated(buf_i8_kind)) deallocate(buf_i8_kind)
301  if (allocated(buf_r4_kind)) deallocate(buf_r4_kind)
302  if (allocated(buf_r8_kind)) deallocate(buf_r8_kind)
303  endif
304 end subroutine compressed_write_2d
305 
306 
307 !> @brief For variables without a compressed dimension, this routine simply wraps
308 !! netcdf_write data. If the variable does have a compressed axis, the I/O root
309 !! gathers the data from the rest of the ranks and then writes the combined data to
310 !! the netcdf file.
311 subroutine compressed_write_3d(fileobj, variable_name, cdata, unlim_dim_level, &
312  corner, edge_lengths)
313 
314  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
315  character(len=*), intent(in) :: variable_name !< Variable name.
316  class(*), contiguous, intent(in), target :: cdata(:,:,:) !< Compressed data that
317  !! will be gathered and
318  !! written to the
319  !! netcdf file.
320  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
321  !! dimension.
322  integer, dimension(3), intent(in), optional :: corner !< Array of starting
323  !! indices describing
324  !! where the data
325  !! will be written to.
326  integer, dimension(3), intent(in), optional :: edge_lengths !< The number of
327  !! elements that
328  !! will be written
329  !! in each dimension.
330 
331  integer, dimension(2) :: compressed_dim_index !< index of the compressed dimension
332  !! compressed_dim_index(1) relative to cdata
333  !! compressed_dim_index(2) relative to the fileobj
334  integer, dimension(3) :: c !! corners of the data to read
335  integer, dimension(3) :: e !< "edges" number of points to read
336 
337  integer(kind=i4_kind), dimension(:,:,:), allocatable :: buf_i4_kind !< Global buffer of data
338  integer(kind=i8_kind), dimension(:,:,:), allocatable :: buf_i8_kind !< Global buffer of data
339  real (kind=r4_kind), dimension(:,:,:), allocatable :: buf_r4_kind !< Global buffer of data
340  real (kind=r8_kind), dimension(:,:,:), allocatable :: buf_r8_kind !< Global buffer of data
341 
342  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
343  class(*), pointer :: cdata_dummy(:,:,:,:)
344 
345  integer :: index(1) !< index of the PE in the pelist
346 
347  integer :: is !< Starting index of the first dimension
348  integer :: ie !< Ending index of the first dimension
349  integer :: js !< Starting index of the second dimension
350  integer :: je !< Ending index of the second dimension
351 
352  append_error_msg = "compressed_write_3d: file:"//trim(fileobj%path)//" variable:"//trim(variable_name)
353 
354  compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
355  if (compressed_dim_index(1) .eq. dimension_not_found) then
356  call netcdf_write_data(fileobj, variable_name, cdata, &
357  unlim_dim_level=unlim_dim_level, corner=corner, &
358  edge_lengths=edge_lengths)
359  return
360  else if (compressed_dim_index(1) .eq. 3) then
361  cdata_dummy(1:size(cdata,1), 1:size(cdata,2), 1:size(cdata,3), 1:1) => cdata(:,:,:)
362  call compressed_write_4d(fileobj, variable_name, cdata_dummy, unlim_dim_level)
363  endif
364 
365  e(:) = shape(cdata)
366  !The root pe creates a buffer big enough to store the data:
367  if (fileobj%is_root) then
368  e(compressed_dim_index(1)) = sum(fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems)
369  select type(cdata)
370  type is (integer(kind=i4_kind))
371  call allocate_array(buf_i4_kind, e)
372  type is (integer(kind=i8_kind))
373  call allocate_array(buf_i8_kind, e)
374  type is (real(kind=r4_kind))
375  call allocate_array(buf_r4_kind, e)
376  type is (real(kind=r8_kind))
377  call allocate_array(buf_r8_kind, e)
378  class default
379  call error("unsupported variable type: "//trim(append_error_msg))
380  end select
381  endif
382 
383  c(:) = 1
384  index = findloc(fileobj%pelist, mpp_pe())
385 
386  c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(index(1))
387  e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(index(1))
388 
389  if (compressed_dim_index(1) .eq. 1) then
390  is = c(compressed_dim_index(1))
391  ie = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
392  js = 1
393  je = size(cdata,2)
394  else
395  js = c(compressed_dim_index(1))
396  je = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
397  is = 1
398  ie = size(cdata,2)
399  endif
400 
401  select type(cdata)
402  type is (integer(kind=i4_kind))
403  call mpp_gather(is, ie, js, je, size(cdata,3), &
404  fileobj%pelist, cdata, buf_i4_kind, fileobj%is_root)
405  if (fileobj%is_root) call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
406  unlim_dim_level=unlim_dim_level)
407  type is (integer(kind=i8_kind))
408  call mpp_gather(is, ie, js, je, size(cdata,3), &
409  fileobj%pelist, cdata, buf_i8_kind, fileobj%is_root)
410  if (fileobj%is_root) call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
411  unlim_dim_level=unlim_dim_level)
412  type is (real(kind=r4_kind))
413  call mpp_gather(is, ie, js, je, size(cdata,3), &
414  fileobj%pelist, cdata, buf_r4_kind, fileobj%is_root)
415  if (fileobj%is_root) call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
416  unlim_dim_level=unlim_dim_level)
417  type is (real(kind=r8_kind))
418  call mpp_gather(is, ie, js, je, size(cdata,3), &
419  fileobj%pelist, cdata, buf_r8_kind, fileobj%is_root)
420  if (fileobj%is_root) call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
421  unlim_dim_level=unlim_dim_level)
422  class default
423  call error("unsupported variable type: "//trim(append_error_msg))
424  end select
425 
426  if (fileobj%is_root) then
427  if (allocated(buf_i4_kind)) deallocate(buf_i4_kind)
428  if (allocated(buf_i8_kind)) deallocate(buf_i8_kind)
429  if (allocated(buf_r4_kind)) deallocate(buf_r4_kind)
430  if (allocated(buf_r8_kind)) deallocate(buf_r8_kind)
431  endif
432 end subroutine compressed_write_3d
433 
434 
435 !> @brief For variables without a compressed dimension, this routine simply wraps
436 !! netcdf_write data. If the variable does have a compressed axis, the I/O root
437 !! gathers the data from the rest of the ranks and then writes the combined data to
438 !! the netcdf file.
439 subroutine compressed_write_4d(fileobj, variable_name, cdata, unlim_dim_level, &
440  corner, edge_lengths)
441 
442  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
443  character(len=*), intent(in) :: variable_name !< Variable name.
444  class(*), dimension(:,:,:,:), intent(in) :: cdata !< Compressed data that
445  !! will be gathered and
446  !! written to the
447  !! netcdf file.
448  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
449  !! dimension.
450  integer, dimension(4), intent(in), optional :: corner !< Array of starting
451  !! indices describing
452  !! where the data
453  !! will be written to.
454  integer, dimension(4), intent(in), optional :: edge_lengths !< The number of
455  !! elements that
456  !! will be written
457  !! in each dimension.
458 
459  integer, dimension(2) :: compressed_dim_index
460  integer, dimension(4) :: c
461  integer, dimension(4) :: e
462  integer :: i
463  integer(kind=i4_kind), dimension(:,:,:,:), allocatable :: buf_i4_kind
464  integer(kind=i8_kind), dimension(:,:,:,:), allocatable :: buf_i8_kind
465  real(kind=r4_kind), dimension(:,:,:,:), allocatable :: buf_r4_kind
466  real(kind=r8_kind), dimension(:,:,:,:), allocatable :: buf_r8_kind
467  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
468 
469  append_error_msg = "compressed_write_4d: file:"//trim(fileobj%path)//" variable:"//trim(variable_name)
470 
471  compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
472  if (compressed_dim_index(1) .eq. dimension_not_found) then
473  call netcdf_write_data(fileobj, variable_name, cdata, &
474  unlim_dim_level=unlim_dim_level, corner=corner, &
475  edge_lengths=edge_lengths)
476  return
477  endif
478 
479  !Gather the data onto the I/O root and write it out.
480  if (fileobj%is_root) then
481  c(:) = 1
482  e(:) = shape(cdata)
483  do i = 1, size(fileobj%pelist)
484  c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(i)
485  e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(i)
486  if (i .eq. 1) then
487  call netcdf_write_data(fileobj, variable_name, cdata, &
488  unlim_dim_level=unlim_dim_level, corner=c, &
489  edge_lengths=e)
490  else
491  select type(cdata)
492  type is (integer(kind=i4_kind))
493  call allocate_array(buf_i4_kind, e)
494  call mpp_recv(buf_i4_kind, size(buf_i4_kind), fileobj%pelist(i), block=.true.)
495  call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
496  unlim_dim_level=unlim_dim_level, corner=c, &
497  edge_lengths=e)
498  deallocate(buf_i4_kind)
499  type is (integer(kind=i8_kind))
500  call allocate_array(buf_i8_kind, e)
501  call mpp_recv(buf_i8_kind, size(buf_i8_kind), fileobj%pelist(i), block=.true.)
502  call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
503  unlim_dim_level=unlim_dim_level, corner=c, &
504  edge_lengths=e)
505  deallocate(buf_i8_kind)
506  type is (real(kind=r4_kind))
507  call allocate_array(buf_r4_kind, e)
508  call mpp_recv(buf_r4_kind, size(buf_r4_kind), fileobj%pelist(i), block=.true.)
509  call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
510  unlim_dim_level=unlim_dim_level, corner=c, &
511  edge_lengths=e)
512  deallocate(buf_r4_kind)
513  type is (real(kind=r8_kind))
514  call allocate_array(buf_r8_kind, e)
515  call mpp_recv(buf_r8_kind, size(buf_r8_kind), fileobj%pelist(i), block=.true.)
516  call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
517  unlim_dim_level=unlim_dim_level, corner=c, &
518  edge_lengths=e)
519  deallocate(buf_r8_kind)
520  class default
521  call error("unsupported variable type: "//trim(append_error_msg))
522  end select
523  endif
524  enddo
525  else
526  select type(cdata)
527  type is (integer(kind=i4_kind))
528  call mpp_send(cdata, size(cdata), fileobj%io_root)
529  type is (integer(kind=i8_kind))
530  call mpp_send(cdata, size(cdata), fileobj%io_root)
531  type is (real(kind=r4_kind))
532  call mpp_send(cdata, size(cdata), fileobj%io_root)
533  type is (real(kind=r8_kind))
534  call mpp_send(cdata, size(cdata), fileobj%io_root)
535  class default
536  call error("unsupported variable type: "//trim(append_error_msg))
537  end select
538  call mpp_sync_self(check=event_send)
539  endif
540 end subroutine compressed_write_4d
541 
542 
543 !> @brief For variables without a compressed dimension, this routine simply wraps
544 !! netcdf_write data. If the variable does have a compressed axis, the I/O root
545 !! gathers the data from the rest of the ranks and then writes the combined data to
546 !! the netcdf file.
547 subroutine compressed_write_5d(fileobj, variable_name, cdata, unlim_dim_level, &
548  corner, edge_lengths)
549 
550  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
551  character(len=*), intent(in) :: variable_name !< Variable name.
552  class(*), dimension(:,:,:,:,:), intent(in) :: cdata !< Compressed data that
553  !! will be gathered and
554  !! written to the
555  !! netcdf file.
556  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
557  !! dimension.
558  integer, dimension(5), intent(in), optional :: corner !< Array of starting
559  !! indices describing
560  !! where the data
561  !! will be written to.
562  integer, dimension(5), intent(in), optional :: edge_lengths !< The number of
563  !! elements that
564  !! will be written
565  !! in each dimension.
566 
567  integer, dimension(2) :: compressed_dim_index
568  integer, dimension(5) :: c
569  integer, dimension(5) :: e
570  integer :: i
571  integer(kind=i4_kind), dimension(:,:,:,:,:), allocatable :: buf_i4_kind
572  integer(kind=i8_kind), dimension(:,:,:,:,:), allocatable :: buf_i8_kind
573  real(kind=r4_kind), dimension(:,:,:,:,:), allocatable :: buf_r4_kind
574  real(kind=r8_kind), dimension(:,:,:,:,:), allocatable :: buf_r8_kind
575  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
576 
577  append_error_msg = "compressed_write_5d: file:"//trim(fileobj%path)//" variable:"//trim(variable_name)
578 
579  compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
580  if (compressed_dim_index(1) .eq. dimension_not_found) then
581  call netcdf_write_data(fileobj, variable_name, cdata, &
582  unlim_dim_level=unlim_dim_level, corner=corner, &
583  edge_lengths=edge_lengths)
584  return
585  endif
586 
587  !Gather the data onto the I/O root and write it out.
588  if (fileobj%is_root) then
589  c(:) = 1
590  e(:) = shape(cdata)
591  do i = 1, size(fileobj%pelist)
592  c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(i)
593  e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(i)
594  if (i .eq. 1) then
595  call netcdf_write_data(fileobj, variable_name, cdata, &
596  unlim_dim_level=unlim_dim_level, corner=c, &
597  edge_lengths=e)
598  else
599  select type(cdata)
600  type is (integer(kind=i4_kind))
601  call allocate_array(buf_i4_kind, e)
602  call mpp_recv(buf_i4_kind, size(buf_i4_kind), fileobj%pelist(i), block=.true.)
603  call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
604  unlim_dim_level=unlim_dim_level, corner=c, &
605  edge_lengths=e)
606  deallocate(buf_i4_kind)
607  type is (integer(kind=i8_kind))
608  call allocate_array(buf_i8_kind, e)
609  call mpp_recv(buf_i8_kind, size(buf_i8_kind), fileobj%pelist(i), block=.true.)
610  call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
611  unlim_dim_level=unlim_dim_level, corner=c, &
612  edge_lengths=e)
613  deallocate(buf_i8_kind)
614  type is (real(kind=r4_kind))
615  call allocate_array(buf_r4_kind, e)
616  call mpp_recv(buf_r4_kind, size(buf_r4_kind), fileobj%pelist(i), block=.true.)
617  call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
618  unlim_dim_level=unlim_dim_level, corner=c, &
619  edge_lengths=e)
620  deallocate(buf_r4_kind)
621  type is (real(kind=r8_kind))
622  call allocate_array(buf_r8_kind, e)
623  call mpp_recv(buf_r8_kind, size(buf_r8_kind), fileobj%pelist(i), block=.true.)
624  call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
625  unlim_dim_level=unlim_dim_level, corner=c, &
626  edge_lengths=e)
627  deallocate(buf_r8_kind)
628  class default
629  call error("unsupported variable type: "//trim(append_error_msg))
630  end select
631  endif
632  enddo
633  else
634  select type(cdata)
635  type is (integer(kind=i4_kind))
636  call mpp_send(cdata, size(cdata), fileobj%io_root)
637  type is (integer(kind=i8_kind))
638  call mpp_send(cdata, size(cdata), fileobj%io_root)
639  type is (real(kind=r4_kind))
640  call mpp_send(cdata, size(cdata), fileobj%io_root)
641  type is (real(kind=r8_kind))
642  call mpp_send(cdata, size(cdata), fileobj%io_root)
643  class default
644  call error("unsupported variable type: "//trim(append_error_msg))
645  end select
646  call mpp_sync_self(check=event_send)
647  endif
648 end subroutine compressed_write_5d
649 
650 
651 !> @brief Wrapper to distinguish interfaces.
652 subroutine compressed_write_1d_wrap(fileobj, variable_name, cdata, unlim_dim_level, &
653  corner, edge_lengths)
654 
655  type(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
656  character(len=*), intent(in) :: variable_name !< Variable name.
657  class(*), dimension(:), intent(in) :: cdata !< Compressed data that
658  !! will be gathered and
659  !! written to the
660  !! netcdf file.
661  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
662  !! dimension.
663  integer, dimension(1), intent(in), optional :: corner !< Array of starting
664  !! indices describing
665  !! where the data
666  !! will be written to.
667  integer, dimension(1), intent(in), optional :: edge_lengths !< The number of
668  !! elements that
669  !! will be written
670  !! in each dimension.
671 
672  call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
673  edge_lengths=edge_lengths)
674 
675 end subroutine compressed_write_1d_wrap
676 
677 
678 !> @brief Wrapper to distinguish interfaces.
679 subroutine compressed_write_2d_wrap(fileobj, variable_name, cdata, unlim_dim_level, &
680  corner, edge_lengths)
681 
682  type(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
683  character(len=*), intent(in) :: variable_name !< Variable name.
684  class(*), dimension(:,:), intent(in) :: cdata !< Compressed data that
685  !! will be gathered and
686  !! written to the
687  !! netcdf file.
688  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
689  !! dimension.
690  integer, dimension(2), intent(in), optional :: corner !< Array of starting
691  !! indices describing
692  !! where the data
693  !! will be written to.
694  integer, dimension(2), intent(in), optional :: edge_lengths !< The number of
695  !! elements that
696  !! will be written
697  !! in each dimension.
698 
699  call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
700  edge_lengths=edge_lengths)
701 end subroutine compressed_write_2d_wrap
702 
703 
704 !> @brief Wrapper to distinguish interfaces.
705 subroutine compressed_write_3d_wrap(fileobj, variable_name, cdata, unlim_dim_level, &
706  corner, edge_lengths)
707 
708  type(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
709  character(len=*), intent(in) :: variable_name !< Variable name.
710  class(*), dimension(:,:,:), intent(in) :: cdata !< Compressed data that
711  !! will be gathered and
712  !! written to the
713  !! netcdf file.
714  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
715  !! dimension.
716  integer, dimension(3), intent(in), optional :: corner !< Array of starting
717  !! indices describing
718  !! where the data
719  !! will be written to.
720  integer, dimension(3), intent(in), optional :: edge_lengths !< The number of
721  !! elements that
722  !! will be written
723  !! in each dimension.
724 
725  call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
726  edge_lengths=edge_lengths)
727 end subroutine compressed_write_3d_wrap
728 
729 
730 !> @brief Wrapper to distinguish interfaces.
731 subroutine compressed_write_4d_wrap(fileobj, variable_name, cdata, unlim_dim_level, &
732  corner, edge_lengths)
733 
734  type(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
735  character(len=*), intent(in) :: variable_name !< Variable name.
736  class(*), dimension(:,:,:,:), intent(in) :: cdata !< Compressed data that
737  !! will be gathered and
738  !! written to the
739  !! netcdf file.
740  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
741  !! dimension.
742  integer, dimension(4), intent(in), optional :: corner !< Array of starting
743  !! indices describing
744  !! where the data
745  !! will be written to.
746  integer, dimension(4), intent(in), optional :: edge_lengths !< The number of
747  !! elements that
748  !! will be written
749  !! in each dimension.
750 
751  call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
752  edge_lengths=edge_lengths)
753 end subroutine compressed_write_4d_wrap
754 
755 
756 !> @brief Wrapper to distinguish interfaces.
757 subroutine compressed_write_5d_wrap(fileobj, variable_name, cdata, unlim_dim_level, &
758  corner, edge_lengths)
759 
760  type(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
761  character(len=*), intent(in) :: variable_name !< Variable name.
762  class(*), dimension(:,:,:,:,:), intent(in) :: cdata !< Compressed data that
763  !! will be gathered and
764  !! written to the
765  !! netcdf file.
766  integer, intent(in), optional :: unlim_dim_level !< Level for the unlimited
767  !! dimension.
768  integer, dimension(5), intent(in), optional :: corner !< Array of starting
769  !! indices describing
770  !! where the data
771  !! will be written to.
772  integer, dimension(5), intent(in), optional :: edge_lengths !< The number of
773  !! elements that
774  !! will be written
775  !! in each dimension.
776 
777  call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
778  edge_lengths=edge_lengths)
779 end subroutine compressed_write_5d_wrap
780 !> @}
subroutine compressed_write_1d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine compressed_write_3d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_5d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine compressed_write_0d(fileobj, variable_name, cdata, unlim_dim_level, corner)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_5d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_0d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner)
Wrapper to distinguish interfaces.
subroutine compressed_write_4d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_2d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine compressed_write_1d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_4d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine compressed_write_2d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_3d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
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...
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407