27 subroutine char_write_0d(fileobj, variable_name, variable_data, append_error_msg, err, varid)
29 class(fmsnetcdffile_t),
intent(in) :: fileobj
30 character(len=*),
intent(in) :: variable_name
31 character(len=*),
intent(in) :: variable_data
32 character(len=200),
intent(in) :: append_error_msg
33 integer,
intent(inout) :: err
34 integer,
intent(in) :: varid
37 integer,
dimension(:),
allocatable :: start
38 integer,
dimension(:),
allocatable :: dimsizes
39 character,
dimension(:),
allocatable :: charbuf
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))
47 allocate(start(ndims))
49 allocate(dimsizes(ndims))
50 call get_variable_size(fileobj, variable_name, dimsizes, .false.)
51 call allocate_array(charbuf, dimsizes)
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))
58 charbuf(i) = variable_data(i:i)
60 err = nf90_put_var(fileobj%ncid, varid, charbuf, start=start, count=dimsizes)
67 subroutine character_write_1d(fileobj, variable_name, variable_data, append_error_msg, err, varid)
69 class(fmsnetcdffile_t),
intent(in) :: fileobj
70 character(len=*),
intent(in) :: variable_name
71 character(len=*),
dimension(:),
intent(in) :: variable_data
72 character(len=200),
intent(in) :: append_error_msg
73 integer,
intent(inout) :: err
74 integer,
intent(in) :: varid
77 integer,
dimension(:),
allocatable :: start
78 integer,
dimension(:),
allocatable :: dimsizes
79 character,
dimension(:,:),
allocatable :: charbuf
80 character(len=1024) :: sbuf
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))
89 allocate(start(ndims))
91 allocate(dimsizes(ndims))
92 call get_variable_size(fileobj, variable_name, dimsizes, .false.)
93 call allocate_array(charbuf, dimsizes)
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))
99 if (
size(variable_data) .ne. dimsizes(2))
then
100 call error(
"incorrect size of variable_data array: "//trim(append_error_msg))
102 do j = 1, dimsizes(2)
103 call string_copy(sbuf, variable_data(j))
105 charbuf(i,j) = sbuf(i:i)
108 err = nf90_put_var(fileobj%ncid, varid, charbuf, start=start, count=dimsizes)
118 class(fmsnetcdffile_t),
intent(in) :: fileobj
119 character(len=*),
intent(in) :: variable_name
120 class(*),
intent(in) :: variable_data
121 integer,
intent(in),
optional :: unlim_dim_level
122 integer,
intent(in),
optional :: corner
129 integer :: unlim_dim_index
130 integer,
dimension(1) :: c
131 character(len=200) :: append_error_msg
133 append_error_msg =
"netcdf_write_data_0d: file:"//trim(fileobj%path)//
" variable: "//trim(variable_name)
134 if (fileobj%is_root)
then
136 if (
present(corner))
then
139 if (
present(unlim_dim_level))
then
140 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
142 if (unlim_dim_index .ne. 1)
then
143 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
145 c(unlim_dim_index) = unlim_dim_level
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)
161 call error(
"Unsupported variable type: "//trim(append_error_msg))
163 call check_netcdf_code(err, append_error_msg)
170 corner, edge_lengths)
172 class(fmsnetcdffile_t),
intent(in) :: fileobj
173 character(len=*),
intent(in) :: variable_name
174 class(*),
dimension(:),
intent(in) :: variable_data
175 integer,
intent(in),
optional :: unlim_dim_level
176 integer,
dimension(1),
intent(in),
optional :: corner
180 integer,
dimension(1),
intent(in),
optional :: edge_lengths
187 integer :: unlim_dim_index
188 integer,
dimension(2) :: c
189 integer,
dimension(2) :: e
190 character(len=200) :: append_error_msg
192 append_error_msg =
"netcdf_write_data_1d: file:"//trim(fileobj%path)//
" variable: "//trim(variable_name)
194 if (fileobj%is_root)
then
196 if (
present(corner))
then
200 if (
present(edge_lengths))
then
201 e(1:1) = edge_lengths(:)
203 e(1:1) = shape(variable_data)
205 if (
present(unlim_dim_level))
then
206 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
208 if (unlim_dim_index .ne. 2)
then
209 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
211 c(unlim_dim_index) = unlim_dim_level
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)
227 call error(
"Unsupported variable type: "//trim(append_error_msg))
229 call check_netcdf_code(err, append_error_msg)
236 corner, edge_lengths)
238 class(fmsnetcdffile_t),
intent(in) :: fileobj
239 character(len=*),
intent(in) :: variable_name
240 class(*),
dimension(:,:),
intent(in) :: variable_data
241 integer,
intent(in),
optional :: unlim_dim_level
242 integer,
dimension(2),
intent(in),
optional :: corner
246 integer,
dimension(2),
intent(in),
optional :: edge_lengths
253 integer :: unlim_dim_index
254 integer,
dimension(3) :: c
255 integer,
dimension(3) :: e
256 character(len=200) :: append_error_msg
258 append_error_msg =
"netcdf_write_data_2d: file:"//trim(fileobj%path)//
" variable: "//trim(variable_name)
260 if (fileobj%is_root)
then
262 if (
present(corner))
then
266 if (
present(edge_lengths))
then
267 e(1:2) = edge_lengths(:)
269 e(1:2) = shape(variable_data)
271 if (
present(unlim_dim_level))
then
272 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
274 if (unlim_dim_index .ne. 3)
then
275 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
277 c(unlim_dim_index) = unlim_dim_level
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)
291 call error(
"Unsupported variable type: "//trim(append_error_msg))
293 call check_netcdf_code(err, append_error_msg)
300 corner, edge_lengths)
302 class(fmsnetcdffile_t),
intent(in) :: fileobj
303 character(len=*),
intent(in) :: variable_name
304 class(*),
dimension(:,:,:),
intent(in) :: variable_data
305 integer,
intent(in),
optional :: unlim_dim_level
307 integer,
dimension(3),
intent(in),
optional :: corner
311 integer,
dimension(3),
intent(in),
optional :: edge_lengths
318 integer :: unlim_dim_index
319 integer,
dimension(4) :: c
320 integer,
dimension(4) :: e
321 character(len=200) :: append_error_msg
323 append_error_msg =
"netcdf_write_data_3d: file:"//trim(fileobj%path)//
" variable: "//trim(variable_name)
325 if (fileobj%is_root)
then
327 if (
present(corner))
then
331 if (
present(edge_lengths))
then
332 e(1:3) = edge_lengths(:)
334 e(1:3) = shape(variable_data)
336 if (
present(unlim_dim_level))
then
337 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
339 if (unlim_dim_index .ne. 4)
then
340 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
342 c(unlim_dim_index) = unlim_dim_level
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)
356 call error(
"Unsupported variable type: "//trim(append_error_msg))
358 call check_netcdf_code(err, append_error_msg)
365 corner, edge_lengths)
367 class(fmsnetcdffile_t),
intent(in) :: fileobj
368 character(len=*),
intent(in) :: variable_name
369 class(*),
dimension(:,:,:,:),
intent(in) :: variable_data
370 integer,
intent(in),
optional :: unlim_dim_level
371 integer,
dimension(4),
intent(in),
optional :: corner
375 integer,
dimension(4),
intent(in),
optional :: edge_lengths
382 integer :: unlim_dim_index
383 integer,
dimension(5) :: c
384 integer,
dimension(5) :: e
385 character(len=200) :: append_error_msg
387 append_error_msg =
"netcdf_write_data_4d: file:"//trim(fileobj%path)//
" variable: "//trim(variable_name)
389 if (fileobj%is_root)
then
391 if (
present(corner))
then
395 if (
present(edge_lengths))
then
396 e(1:4) = edge_lengths(:)
398 e(1:4) = shape(variable_data)
400 if (
present(unlim_dim_level))
then
401 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
403 if (unlim_dim_index /= -1)
then
404 c(unlim_dim_index) = unlim_dim_level
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)
419 call error(
"Unsupported variable type: "//trim(append_error_msg))
421 call check_netcdf_code(err, append_error_msg)
428 corner, edge_lengths)
430 class(fmsnetcdffile_t),
intent(in) :: fileobj
431 character(len=*),
intent(in) :: variable_name
432 class(*),
dimension(:,:,:,:,:),
intent(in),
target :: variable_data
434 integer,
intent(in),
optional :: unlim_dim_level
436 integer,
dimension(5),
intent(in),
optional :: corner
440 integer,
dimension(5),
intent(in),
optional :: edge_lengths
447 integer :: unlim_dim_index
448 integer,
dimension(6) :: c
449 integer,
dimension(6) :: e
450 character(len=200) :: append_error_msg
452 append_error_msg =
"netcdf_write_data_5d: file:"//trim(fileobj%path)//
" variable: "//trim(variable_name)
454 if (fileobj%is_root)
then
456 if (
present(corner))
then
460 if (
present(edge_lengths))
then
461 e(1:5) = edge_lengths(:)
463 e(1:5) = shape(variable_data)
465 if (
present(unlim_dim_level))
then
466 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
468 if (unlim_dim_index .ne. 6)
then
469 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
471 c(unlim_dim_index) = unlim_dim_level
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)
485 call error(
"Unsupported variable type: "//trim(append_error_msg))
487 call check_netcdf_code(err, append_error_msg)
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.