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