FMS  2025.04
Flexible Modeling System
netcdf_write_data.inc
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 !> @file
19 !> @brief Routines for writing variable values to netcdf files with different
20 !! dimension sizes for the @ref netcdf_write_data interface
21 
22 !> @addtogroup netcdf_io_mod
23 !> @{
24 
25 !> @brief Character write 0d function.
26 subroutine char_write_0d(fileobj, variable_name, variable_data, append_error_msg, err, varid)
27 
28  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
29  character(len=*), intent(in) :: variable_name !< Variable name.
30  character(len=*), intent(in) :: variable_data !< Data that will be written.
31  character(len=200), intent(in) :: append_error_msg !< Msg to be appended to FATAL error message
32  integer, intent(inout) :: err
33  integer, intent(in) :: varid
34 
35  integer :: ndims
36  integer, dimension(:), allocatable :: start
37  integer, dimension(:), allocatable :: dimsizes
38  character, dimension(:), allocatable :: charbuf
39  integer :: i
40  integer :: tlen
41 
42  ndims = get_variable_num_dimensions(fileobj, variable_name, broadcast=.false.)
43  if (ndims .ne. 1) then
44  call error("currently only scalar and 1d character writes are supported: "//trim(append_error_msg))
45  endif
46  allocate(start(ndims))
47  start(:) = 1
48  allocate(dimsizes(ndims))
49  call get_variable_size(fileobj, variable_name, dimsizes, .false.)
50  call allocate_array(charbuf, dimsizes)
51  charbuf(:) = ""
52  tlen = len_trim(variable_data)
53  if (tlen .gt. dimsizes(1)) then
54  call error("character buffer is too big; decrease length: "//trim(append_error_msg))
55  endif
56  do i = 1, tlen
57  charbuf(i) = variable_data(i:i)
58  enddo
59  err = nf90_put_var(fileobj%ncid, varid, charbuf, start=start, count=dimsizes)
60  deallocate(charbuf)
61  deallocate(dimsizes)
62  deallocate(start)
63 end subroutine char_write_0d
64 
65 !> @brief Character write 1d function.
66 subroutine character_write_1d(fileobj, variable_name, variable_data, append_error_msg, err, varid)
67 
68  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
69  character(len=*), intent(in) :: variable_name !< Variable name.
70  character(len=*), dimension(:), intent(in) :: variable_data !< Data that will be written.
71  character(len=200), intent(in) :: append_error_msg !< Msg to be appended to FATAL error message
72  integer, intent(inout) :: err
73  integer, intent(in) :: varid
74 
75  integer :: ndims
76  integer, dimension(:), allocatable :: start
77  integer, dimension(:), allocatable :: dimsizes
78  character, dimension(:,:), allocatable :: charbuf
79  character(len=1024) :: sbuf
80  integer :: i
81  integer :: j
82  integer :: tlen
83 
84  ndims = get_variable_num_dimensions(fileobj, variable_name, broadcast=.false.)
85  if (ndims .ne. 2) then
86  call error("currently only scalar and 1d character writes are supported: "//trim(append_error_msg))
87  endif
88  allocate(start(ndims))
89  start(:) = 1
90  allocate(dimsizes(ndims))
91  call get_variable_size(fileobj, variable_name, dimsizes, .false.)
92  call allocate_array(charbuf, dimsizes)
93  charbuf(:,:) = ""
94  tlen = len(variable_data(1))
95  if (tlen .gt. dimsizes(1)) then
96  call error("character buffer is too big; decrease length: "//trim(append_error_msg))
97  endif
98  if (size(variable_data) .ne. dimsizes(2)) then
99  call error("incorrect size of variable_data array: "//trim(append_error_msg))
100  endif
101  do j = 1, dimsizes(2)
102  call string_copy(sbuf, variable_data(j))
103  do i = 1, tlen
104  charbuf(i,j) = sbuf(i:i)
105  enddo
106  enddo
107  err = nf90_put_var(fileobj%ncid, varid, charbuf, start=start, count=dimsizes)
108  deallocate(charbuf)
109  deallocate(dimsizes)
110  deallocate(start)
111 end subroutine character_write_1d
112 
113 !> @brief Write data to a variable in a netcdf file.
114 subroutine netcdf_write_data_0d(fileobj, variable_name, variable_data, unlim_dim_level, &
115  corner)
116 
117  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
118  character(len=*), intent(in) :: variable_name !< Variable name.
119  class(*), intent(in) :: variable_data !< Data that will be written.
120  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension level.
121  integer, intent(in), optional :: corner !< Array of starting
122  !! indices describing
123  !! where the data
124  !! will be written to.
125 
126  integer :: err
127  integer :: varid
128  integer :: unlim_dim_index
129  integer, dimension(1) :: c
130  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
131 
132  append_error_msg = "netcdf_write_data_0d: file:"//trim(fileobj%path)//" variable: "//trim(variable_name)
133  if (fileobj%is_root) then
134  c(:) = 1
135  if (present(corner)) then
136  c(1) = corner
137  endif
138  if (present(unlim_dim_level)) then
139  unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
140  broadcast=.false.)
141  if (unlim_dim_index .ne. 1) then
142  call error("unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
143  endif
144  c(unlim_dim_index) = unlim_dim_level
145  endif
146  call set_netcdf_mode(fileobj%ncid, data_mode)
147  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
148  select type(variable_data)
149  type is (integer(kind=i4_kind))
150  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c)
151  type is (integer(kind=i8_kind))
152  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c)
153  type is (real(kind=r4_kind))
154  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c)
155  type is (real(kind=r8_kind))
156  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c)
157  type is (character(len=*))
158  call char_write_0d(fileobj, variable_name, variable_data, append_error_msg, err, varid)
159  class default
160  call error("Unsupported variable type: "//trim(append_error_msg))
161  end select
162  call check_netcdf_code(err, append_error_msg)
163  endif
164 end subroutine netcdf_write_data_0d
165 
166 
167 !> @brief Write data to a variable in a netcdf file.
168 subroutine netcdf_write_data_1d(fileobj, variable_name, variable_data, unlim_dim_level, &
169  corner, edge_lengths)
170 
171  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
172  character(len=*), intent(in) :: variable_name !< Variable name.
173  class(*), dimension(:), intent(in) :: variable_data !< Data that will be written.
174  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension level.
175  integer, dimension(1), intent(in), optional :: corner !< Array of starting
176  !! indices describing
177  !! where the data
178  !! will be written to.
179  integer, dimension(1), intent(in), optional :: edge_lengths !< The number of
180  !! elements that
181  !! will be written
182  !! in each dimension.
183 
184  integer :: err
185  integer :: varid
186  integer :: unlim_dim_index
187  integer, dimension(2) :: c
188  integer, dimension(2) :: e
189  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
190 
191  append_error_msg = "netcdf_write_data_1d: file:"//trim(fileobj%path)//" variable: "//trim(variable_name)
192 
193  if (fileobj%is_root) then
194  c(:) = 1
195  if (present(corner)) then
196  c(1:1) = corner(:)
197  endif
198  e(:) = 1
199  if (present(edge_lengths)) then
200  e(1:1) = edge_lengths(:)
201  else
202  e(1:1) = shape(variable_data)
203  endif
204  if (present(unlim_dim_level)) then
205  unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
206  broadcast=.false.)
207  if (unlim_dim_index .ne. 2) then
208  call error("unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
209  endif
210  c(unlim_dim_index) = unlim_dim_level
211  endif
212  call set_netcdf_mode(fileobj%ncid, data_mode)
213  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
214  select type(variable_data)
215  type is (integer(kind=i4_kind))
216  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
217  type is (integer(kind=i8_kind))
218  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
219  type is (real(kind=r4_kind))
220  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
221  type is (real(kind=r8_kind))
222  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
223  type is (character(len=*))
224  call character_write_1d(fileobj, variable_name, variable_data, append_error_msg, err, varid)
225  class default
226  call error("Unsupported variable type: "//trim(append_error_msg))
227  end select
228  call check_netcdf_code(err, append_error_msg)
229  endif
230 end subroutine netcdf_write_data_1d
231 
232 
233 !> @brief Write data to a variable in a netcdf file.
234 subroutine netcdf_write_data_2d(fileobj, variable_name, variable_data, unlim_dim_level, &
235  corner, edge_lengths)
236 
237  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
238  character(len=*), intent(in) :: variable_name !< Variable name.
239  class(*), dimension(:,:), intent(in) :: variable_data !< Data that will be written.
240  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension level.
241  integer, dimension(2), intent(in), optional :: corner !< Array of starting
242  !! indices describing
243  !! where the data
244  !! will be written to.
245  integer, dimension(2), intent(in), optional :: edge_lengths !< The number of
246  !! elements that
247  !! will be written
248  !! in each dimension.
249 
250  integer :: err
251  integer :: varid
252  integer :: unlim_dim_index
253  integer,dimension(3) :: c
254  integer, dimension(3) :: e
255  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
256 
257  append_error_msg = "netcdf_write_data_2d: file:"//trim(fileobj%path)//" variable: "//trim(variable_name)
258 
259  if (fileobj%is_root) then
260  c(:) = 1
261  if (present(corner)) then
262  c(1:2) = corner(:)
263  endif
264  e(:) = 1
265  if (present(edge_lengths)) then
266  e(1:2) = edge_lengths(:)
267  else
268  e(1:2) = shape(variable_data)
269  endif
270  if (present(unlim_dim_level)) then
271  unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
272  broadcast=.false.)
273  if (unlim_dim_index .ne. 3) then
274  call error("unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
275  endif
276  c(unlim_dim_index) = unlim_dim_level
277  endif
278  call set_netcdf_mode(fileobj%ncid, data_mode)
279  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
280  select type(variable_data)
281  type is (integer(kind=i4_kind))
282  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
283  type is (integer(kind=i8_kind))
284  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
285  type is (real(kind=r4_kind))
286  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
287  type is (real(kind=r8_kind))
288  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c ,count=e)
289  class default
290  call error("Unsupported variable type: "//trim(append_error_msg))
291  end select
292  call check_netcdf_code(err, append_error_msg)
293  endif
294 end subroutine netcdf_write_data_2d
295 
296 
297 !> @brief Write data to a variable in a netcdf file.
298 subroutine netcdf_write_data_3d(fileobj, variable_name, variable_data, unlim_dim_level, &
299  corner, edge_lengths)
300 
301  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
302  character(len=*), intent(in) :: variable_name !< Variable name.
303  class(*), dimension(:,:,:), intent(in) :: variable_data !< Data that will be written.
304  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension
305  !! level.
306  integer, dimension(3), intent(in), optional :: corner !< Array of starting
307  !! indices describing
308  !! where the data
309  !! will be written to.
310  integer, dimension(3), intent(in), optional :: edge_lengths !< The number of
311  !! elements that
312  !! will be written
313  !! in each dimension.
314 
315  integer :: err
316  integer :: varid
317  integer :: unlim_dim_index
318  integer, dimension(4) :: c
319  integer, dimension(4) :: e
320  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
321 
322  append_error_msg = "netcdf_write_data_3d: file:"//trim(fileobj%path)//" variable: "//trim(variable_name)
323 
324  if (fileobj%is_root) then
325  c(:) = 1
326  if (present(corner)) then
327  c(1:3) = corner(:)
328  endif
329  e(:) = 1
330  if (present(edge_lengths)) then
331  e(1:3) = edge_lengths(:)
332  else
333  e(1:3) = shape(variable_data)
334  endif
335  if (present(unlim_dim_level)) then
336  unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
337  broadcast=.false.)
338  if (unlim_dim_index .ne. 4) then
339  call error("unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
340  endif
341  c(unlim_dim_index) = unlim_dim_level
342  endif
343  call set_netcdf_mode(fileobj%ncid, data_mode)
344  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
345  select type(variable_data)
346  type is (integer(kind=i4_kind))
347  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
348  type is (integer(kind=i8_kind))
349  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
350  type is (real(kind=r4_kind))
351  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
352  type is (real(kind=r8_kind))
353  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
354  class default
355  call error("Unsupported variable type: "//trim(append_error_msg))
356  end select
357  call check_netcdf_code(err, append_error_msg)
358  endif
359 end subroutine netcdf_write_data_3d
360 
361 
362 !> @brief Write data to a variable in a netcdf file.
363 subroutine netcdf_write_data_4d(fileobj, variable_name, variable_data, unlim_dim_level, &
364  corner, edge_lengths)
365 
366  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
367  character(len=*), intent(in) :: variable_name !< Variable name.
368  class(*), dimension(:,:,:,:), intent(in) :: variable_data !< Data that will be written.
369  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension level.
370  integer, dimension(4), intent(in), optional :: corner !< Array of starting
371  !! indices describing
372  !! where the data
373  !! will be written to.
374  integer, dimension(4), intent(in), optional :: edge_lengths !< The number of
375  !! elements that
376  !! will be written
377  !! in each dimension.
378 
379  integer :: err
380  integer :: varid
381  integer :: unlim_dim_index
382  integer,dimension(5) :: c
383  integer, dimension(5) :: e
384  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
385 
386  append_error_msg = "netcdf_write_data_4d: file:"//trim(fileobj%path)//" variable: "//trim(variable_name)
387 
388  if (fileobj%is_root) then
389  c(:) = 1
390  if (present(corner)) then
391  c(1:4) = corner(:)
392  endif
393  e(:) = 1
394  if (present(edge_lengths)) then
395  e(1:4) = edge_lengths(:)
396  else
397  e(1:4) = shape(variable_data)
398  endif
399  if (present(unlim_dim_level)) then
400  unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
401  broadcast=.false.)
402  if (unlim_dim_index /= -1) then
403  c(unlim_dim_index) = unlim_dim_level
404  endif
405  endif
406  call set_netcdf_mode(fileobj%ncid, data_mode)
407  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
408  select type(variable_data)
409  type is (integer(kind=i4_kind))
410  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
411  type is (integer(kind=i8_kind))
412  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
413  type is (real(kind=r4_kind))
414  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
415  type is (real(kind=r8_kind))
416  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
417  class default
418  call error("Unsupported variable type: "//trim(append_error_msg))
419  end select
420  call check_netcdf_code(err, append_error_msg)
421  endif
422 end subroutine netcdf_write_data_4d
423 
424 
425 !> @brief Write data to a variable in a netcdf file.
426 subroutine netcdf_write_data_5d(fileobj, variable_name, variable_data, unlim_dim_level, &
427  corner, edge_lengths)
428 
429  class(fmsnetcdffile_t), intent(in) :: fileobj !< File object.
430  character(len=*), intent(in) :: variable_name !< Variable name.
431  class(*), dimension(:,:,:,:,:), intent(in), target :: variable_data !< Data that will be
432  !! written.
433  integer, intent(in), optional :: unlim_dim_level !< Unlimited dimension
434  !! level.
435  integer, dimension(5), intent(in), optional :: corner !< Array of starting
436  !! indices describing
437  !! where the data
438  !! will be written to.
439  integer, dimension(5), intent(in), optional :: edge_lengths !< The number of
440  !! elements that
441  !! will be written
442  !! in each dimension.
443 
444  integer :: err
445  integer :: varid
446  integer :: unlim_dim_index
447  integer,dimension(6) :: c
448  integer, dimension(6) :: e
449  character(len=200) :: append_error_msg !< Msg to be appended to FATAL error message
450 
451  append_error_msg = "netcdf_write_data_5d: file:"//trim(fileobj%path)//" variable: "//trim(variable_name)
452 
453  if (fileobj%is_root) then
454  c(:) = 1
455  if (present(corner)) then
456  c(1:5) = corner(:)
457  endif
458  e(:) = 1
459  if (present(edge_lengths)) then
460  e(1:5) = edge_lengths(:)
461  else
462  e(1:5) = shape(variable_data)
463  endif
464  if (present(unlim_dim_level)) then
465  unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
466  broadcast=.false.)
467  if (unlim_dim_index .ne. 6) then
468  call error("unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
469  endif
470  c(unlim_dim_index) = unlim_dim_level
471  endif
472  call set_netcdf_mode(fileobj%ncid, data_mode)
473  varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
474  select type(variable_data)
475  type is (integer(kind=i4_kind))
476  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
477  type is (integer(kind=i8_kind))
478  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
479  type is (real(kind=r4_kind))
480  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
481  type is (real(kind=r8_kind))
482  err = nf90_put_var(fileobj%ncid, varid, variable_data, start=c, count=e)
483  class default
484  call error("Unsupported variable type: "//trim(append_error_msg))
485  end select
486  call check_netcdf_code(err, append_error_msg)
487  endif
488 end subroutine netcdf_write_data_5d
subroutine char_write_0d(fileobj, variable_name, variable_data, append_error_msg, err, varid)
Character write 0d function.
subroutine netcdf_write_data_4d(fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
Write data to a variable in a netcdf file.
subroutine character_write_1d(fileobj, variable_name, variable_data, append_error_msg, err, varid)
Character write 1d function.
subroutine netcdf_write_data_0d(fileobj, variable_name, variable_data, unlim_dim_level, corner)
Write data to a variable in a netcdf file.
subroutine netcdf_write_data_3d(fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
Write data to a variable in a netcdf file.
subroutine netcdf_write_data_1d(fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
Write data to a variable in a netcdf file.
subroutine netcdf_write_data_2d(fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
Write data to a variable in a netcdf file.
subroutine netcdf_write_data_5d(fileobj, variable_name, variable_data, unlim_dim_level, corner, edge_lengths)
Write data to a variable in a netcdf file.