26 subroutine char_read_0d(fileobj, variable_name, buf, corner, append_error_msg, err, varid)
27 class(fmsnetcdffile_t),
intent(in) :: fileobj
28 character(len=*),
intent(in) :: variable_name
29 character(len=*),
intent(inout) :: buf
30 integer,
intent(in),
optional :: corner
32 character(len=200),
intent(in):: append_error_msg
33 integer,
intent(inout) :: err
34 integer,
intent(in) :: varid
36 integer,
dimension(2) :: start
38 character,
dimension(:),
allocatable :: charbuf
39 integer,
dimension(:),
allocatable :: dimsizes
43 ndims = get_variable_num_dimensions(fileobj, variable_name, broadcast=.false.)
44 allocate(dimsizes(ndims))
45 call get_variable_size(fileobj, variable_name, dimsizes, broadcast=.false.)
46 allocate(charbuf(dimsizes(1)))
48 if (ndims .eq. 2)
then
49 if (
present(corner))
then
53 elseif (ndims .gt. 2)
then
54 call error(
"Only scalar and 1d string values are currently supported: "//trim(append_error_msg))
56 err = nf90_get_var(fileobj%ncid, varid, charbuf, start=start, count=dimsizes)
57 if (len(buf) .lt. dimsizes(1))
then
58 call error(
"character buffer is too small; increase length: "//trim(append_error_msg))
62 if (charbuf(i) .eq. char(0))
then
72 subroutine char_read_1d(fileobj, variable_name, buf, c, append_error_msg, err, varid)
74 class(fmsnetcdffile_t),
intent(in) :: fileobj
75 character(len=*),
intent(in) :: variable_name
76 character(len=*),
dimension(:),
intent(inout) :: buf
78 integer,
dimension(2),
intent(in) :: c
79 character(len=200),
intent(in) :: append_error_msg
80 integer,
intent(inout) :: err
81 integer,
intent(in) :: varid
84 integer,
dimension(2) :: start
85 integer,
dimension(2) :: dimsizes
86 character,
dimension(:,:),
allocatable :: charbuf
87 character(len=1024) :: sbuf
91 ndims = get_variable_num_dimensions(fileobj, variable_name, broadcast=.false.)
92 if (ndims .ne. 2)
then
93 call error(trim(variable_name)//
"must be 2d dimensional in netcdf file.")
97 call get_variable_size(fileobj, variable_name, dimsizes, .false.)
98 dimsizes(2) = dimsizes(2) - start(2) + 1
99 call allocate_array(charbuf, dimsizes, initialize=.true.)
100 err = nf90_get_var(fileobj%ncid, varid, charbuf, start=start, count=dimsizes)
101 if (len(buf(1)) .lt. dimsizes(1))
then
102 call error(
"character buffer is too small; increase length: "//trim(append_error_msg))
104 if (
size(buf) .lt. dimsizes(2))
then
105 call error(
"incorrect buffer array size:: "//trim(append_error_msg))
107 do i = start(2), start(2)+dimsizes(2)-1
109 do j = 1, dimsizes(1)
110 if (charbuf(j, i-start(2)+1) .eq. char(0))
then
113 sbuf(j:j) = charbuf(j, i-start(2)+1)
115 call string_copy(buf(i-start(2)+1), sbuf)
125 class(fmsnetcdffile_t),
intent(in) :: fileobj
126 character(len=*),
intent(in) :: variable_name
127 class(*),
intent(inout) :: buf
128 integer,
intent(in),
optional :: unlim_dim_level
129 integer,
intent(in),
optional :: corner
131 logical,
intent(in),
optional :: broadcast
141 integer :: unlim_dim_index
142 integer,
dimension(1) :: c
143 character(len=1024),
dimension(1) :: buf1d
144 character(len=200) :: append_error_msg
146 append_error_msg =
"netcdf_read_data_0d: file:"//trim(fileobj%path)//
"- variable:"//trim(variable_name)
148 if (
present(broadcast))
then
153 if (fileobj%use_netcdf_mpi) bcast = .false.
156 if (
present(unlim_dim_level))
then
157 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
159 if (unlim_dim_index .ne. 1)
then
160 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
162 c(unlim_dim_index) = unlim_dim_level
165 if (fileobj%is_root)
then
166 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
168 if (fileobj%use_netcdf_mpi.and.fileobj%use_collective)
then
169 err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
170 call check_netcdf_code(err, append_error_msg)
174 type is (
integer(kind=i4_kind))
175 err = nf90_get_var(fileobj%ncid, varid, buf, start=c)
176 type is (
integer(kind=i8_kind))
177 err = nf90_get_var(fileobj%ncid, varid, buf, start=c)
178 type is (real(kind=r4_kind))
179 err = nf90_get_var(fileobj%ncid, varid, buf, start=c)
180 type is (real(kind=r8_kind))
181 err = nf90_get_var(fileobj%ncid, varid, buf, start=c)
182 type is (
character(len=*))
183 call char_read_0d(fileobj, variable_name, buf, corner, append_error_msg, err, varid)
185 call error(
"Unsupported variable type: "//trim(append_error_msg))
187 call check_netcdf_code(err, append_error_msg)
193 type is (
integer(kind=i4_kind))
194 call mpp_broadcast(buf, fileobj%io_root, pelist=fileobj%pelist)
195 type is (
integer(kind=i8_kind))
196 call mpp_broadcast(buf, fileobj%io_root, pelist=fileobj%pelist)
197 type is (real(kind=r4_kind))
198 call mpp_broadcast(buf, fileobj%io_root, pelist=fileobj%pelist)
199 type is (real(kind=r8_kind))
200 call mpp_broadcast(buf, fileobj%io_root, pelist=fileobj%pelist)
201 type is (
character(len=*))
202 call string_copy(buf1d(1), buf)
203 call mpp_broadcast(buf1d, len(buf1d(1)), fileobj%io_root, pelist=fileobj%pelist)
204 call string_copy(buf, buf1d(1))
206 call error(
"Unsupported variable type: "//trim(append_error_msg))
214 corner, edge_lengths, broadcast)
216 class(fmsnetcdffile_t),
intent(in) :: fileobj
217 character(len=*),
intent(in) :: variable_name
218 class(*),
dimension(:),
intent(inout) :: buf
220 integer,
intent(in),
optional :: unlim_dim_level
222 integer,
dimension(:),
intent(in),
optional :: corner
226 integer,
dimension(:),
intent(in),
optional :: edge_lengths
230 logical,
intent(in),
optional :: broadcast
240 integer :: unlim_dim_index
241 integer,
dimension(2) :: c
242 integer,
dimension(2) :: e
243 character(len=200) :: append_error_msg
245 append_error_msg =
"netcdf_read_data_1d: file:"//trim(fileobj%path)//
"- variable:"//trim(variable_name)
247 if (
present(broadcast))
then
252 if (fileobj%use_netcdf_mpi) bcast = .false.
255 if (
present(corner))
then
259 if (
present(edge_lengths))
then
260 e(1) = edge_lengths(1)
264 if (
present(unlim_dim_level))
then
265 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
267 if (unlim_dim_index .ne. 2)
then
268 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
270 c(unlim_dim_index) = unlim_dim_level
273 if (fileobj%is_root)
then
274 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
276 if (fileobj%use_netcdf_mpi.and.fileobj%use_collective)
then
277 err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
278 call check_netcdf_code(err, append_error_msg)
282 type is (
integer(kind=i4_kind))
283 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
284 type is (
integer(kind=i8_kind))
285 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
286 type is (real(kind=r4_kind))
287 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
288 type is (real(kind=r8_kind))
289 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
290 type is (
character(len=*))
291 call char_read_1d(fileobj, variable_name, buf, c, append_error_msg, err, varid)
293 call error(
"Unsupported variable type: "//trim(append_error_msg))
295 call check_netcdf_code(err, append_error_msg)
301 type is (
integer(kind=i4_kind))
302 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
303 type is (
integer(kind=i8_kind))
304 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
305 type is (real(kind=r4_kind))
306 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
307 type is (real(kind=r8_kind))
308 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
309 type is (
character(len=*))
310 call mpp_broadcast(buf, len(buf(1)), fileobj%io_root, pelist=fileobj%pelist)
312 call error(
"Unsupported variable type: "//trim(append_error_msg))
320 corner, edge_lengths, broadcast)
322 class(fmsnetcdffile_t),
intent(in) :: fileobj
323 character(len=*),
intent(in) :: variable_name
324 class(*),
dimension(:,:),
intent(inout) :: buf
326 integer,
intent(in),
optional :: unlim_dim_level
328 integer,
dimension(:),
intent(in),
optional :: corner
332 integer,
dimension(:),
intent(in),
optional :: edge_lengths
336 logical,
intent(in),
optional :: broadcast
346 integer :: unlim_dim_index
347 integer,
dimension(3) :: c
348 integer,
dimension(3) :: e
349 character(len=200) :: append_error_msg
351 append_error_msg =
"netcdf_read_data_2d: file:"//trim(fileobj%path)//
"- variable:"//trim(variable_name)
353 if (
present(broadcast))
then
358 if (fileobj%use_netcdf_mpi) bcast = .false.
361 if (
present(corner))
then
365 if (
present(edge_lengths))
then
366 e(1:2) = edge_lengths(1:2)
370 if (
present(unlim_dim_level))
then
371 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
373 if (unlim_dim_index .ne. 3)
then
374 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
376 c(unlim_dim_index) = unlim_dim_level
379 if (fileobj%is_root)
then
380 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
382 if (fileobj%use_netcdf_mpi.and.fileobj%use_collective)
then
383 err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
384 call check_netcdf_code(err, append_error_msg)
388 type is (
integer(kind=i4_kind))
389 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
390 type is (
integer(kind=i8_kind))
391 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
392 type is (real(kind=r4_kind))
393 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
394 type is (real(kind=r8_kind))
395 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
397 call error(
"Unsupported variable type: "//trim(append_error_msg))
399 call check_netcdf_code(err, append_error_msg)
405 type is (
integer(kind=i4_kind))
406 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
407 type is (
integer(kind=i8_kind))
408 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
409 type is (real(kind=r4_kind))
410 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
411 type is (real(kind=r8_kind))
412 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
414 call error(
"Unsupported variable type: "//trim(append_error_msg))
422 corner, edge_lengths, broadcast)
424 class(fmsnetcdffile_t),
intent(in) :: fileobj
425 character(len=*),
intent(in) :: variable_name
426 class(*),
dimension(:,:,:),
intent(inout) :: buf
428 integer,
intent(in),
optional :: unlim_dim_level
430 integer,
dimension(:),
intent(in),
optional :: corner
434 integer,
dimension(:),
intent(in),
optional :: edge_lengths
438 logical,
intent(in),
optional :: broadcast
448 integer :: unlim_dim_index
449 integer,
dimension(4) :: c
450 integer,
dimension(4) :: e
451 character(len=200) :: append_error_msg
453 append_error_msg =
"netcdf_read_data_3d: file:"//trim(fileobj%path)//
"- variable:"//trim(variable_name)
455 if (
present(broadcast))
then
460 if (fileobj%use_netcdf_mpi) bcast = .false.
463 if (
present(corner))
then
467 if (
present(edge_lengths))
then
468 e(1:3) = edge_lengths(1:3)
472 if (
present(unlim_dim_level))
then
473 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
475 if (unlim_dim_index .ne. 4)
then
476 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
478 c(unlim_dim_index) = unlim_dim_level
481 if (fileobj%is_root)
then
482 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
484 if (fileobj%use_netcdf_mpi.and.fileobj%use_collective)
then
485 err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
486 call check_netcdf_code(err, append_error_msg)
490 type is (
integer(kind=i4_kind))
491 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
492 type is (
integer(kind=i8_kind))
493 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
494 type is (real(kind=r4_kind))
495 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
496 type is (real(kind=r8_kind))
497 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
499 call error(
"Unsupported variable type: "//trim(append_error_msg))
501 call check_netcdf_code(err, append_error_msg)
507 type is (
integer(kind=i4_kind))
508 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
509 type is (
integer(kind=i8_kind))
510 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
511 type is (real(kind=r4_kind))
512 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
513 type is (real(kind=r8_kind))
514 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
516 call error(
"Unsupported variable type: "//trim(append_error_msg))
524 corner, edge_lengths, broadcast)
526 class(fmsnetcdffile_t),
intent(in) :: fileobj
527 character(len=*),
intent(in) :: variable_name
528 class(*),
dimension(:,:,:,:),
intent(inout) :: buf
530 integer,
intent(in),
optional :: unlim_dim_level
532 integer,
dimension(:),
intent(in),
optional :: corner
536 integer,
dimension(:),
intent(in),
optional :: edge_lengths
540 logical,
intent(in),
optional :: broadcast
550 integer :: unlim_dim_index
551 integer,
dimension(5) :: c
552 integer,
dimension(5) :: e
553 character(len=200) :: append_error_msg
555 append_error_msg =
"netcdf_read_data_4d: file:"//trim(fileobj%path)//
"- variable:"//trim(variable_name)
557 if (
present(broadcast))
then
562 if (fileobj%use_netcdf_mpi) bcast = .false.
565 if (
present(corner))
then
569 if (
present(edge_lengths))
then
570 e(1:4) = edge_lengths(1:4)
574 if (
present(unlim_dim_level))
then
575 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
577 if (unlim_dim_index .ne. 5)
then
578 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
580 c(unlim_dim_index) = unlim_dim_level
583 if (fileobj%is_root)
then
584 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
586 if (fileobj%use_netcdf_mpi.and.fileobj%use_collective)
then
587 err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
588 call check_netcdf_code(err, append_error_msg)
592 type is (
integer(kind=i4_kind))
593 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
594 type is (
integer(kind=i8_kind))
595 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
596 type is (real(kind=r4_kind))
597 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
598 type is (real(kind=r8_kind))
599 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
601 call error(
"Unsupported variable type: "//trim(append_error_msg))
603 call check_netcdf_code(err, append_error_msg)
609 type is (
integer(kind=i4_kind))
610 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
611 type is (
integer(kind=i8_kind))
612 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
613 type is (real(kind=r4_kind))
614 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
615 type is (real(kind=r8_kind))
616 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
618 call error(
"Unsupported variable type: "//trim(append_error_msg))
626 corner, edge_lengths, broadcast)
628 class(fmsnetcdffile_t),
intent(in) :: fileobj
629 character(len=*),
intent(in) :: variable_name
630 class(*),
dimension(:,:,:,:,:),
intent(inout) :: buf
632 integer,
intent(in),
optional :: unlim_dim_level
634 integer,
dimension(:),
intent(in),
optional :: corner
638 integer,
dimension(:),
intent(in),
optional :: edge_lengths
642 logical,
intent(in),
optional :: broadcast
652 integer :: unlim_dim_index
653 integer,
dimension(6) :: c
654 integer,
dimension(6) :: e
655 character(len=200) :: append_error_msg
657 append_error_msg =
"netcdf_read_data_5d: file:"//trim(fileobj%path)//
"- variable:"//trim(variable_name)
659 if (
present(broadcast))
then
664 if (fileobj%use_netcdf_mpi) bcast = .false.
667 if (
present(corner))
then
671 if (
present(edge_lengths))
then
672 e(1:5) = edge_lengths(1:5)
676 if (
present(unlim_dim_level))
then
677 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
679 if (unlim_dim_index .ne. 6)
then
680 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
682 c(unlim_dim_index) = unlim_dim_level
685 if (fileobj%is_root)
then
686 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
688 if (fileobj%use_netcdf_mpi.and.fileobj%use_collective)
then
689 err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
690 call check_netcdf_code(err, append_error_msg)
694 type is (
integer(kind=i4_kind))
695 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
696 type is (
integer(kind=i8_kind))
697 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
698 type is (real(kind=r4_kind))
699 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
700 type is (real(kind=r8_kind))
701 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
703 call error(
"Unsupported variable type: "//trim(append_error_msg))
705 call check_netcdf_code(err, append_error_msg)
711 type is (
integer(kind=i4_kind))
712 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
713 type is (
integer(kind=i8_kind))
714 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
715 type is (real(kind=r4_kind))
716 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
717 type is (real(kind=r8_kind))
718 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
720 call error(
"Unsupported variable type: "//trim(append_error_msg))
subroutine unpack_data_2d(fileobj, varid, varname, var_data)
subroutine unpack_data_1d(fileobj, varid, varname, var_data)
subroutine netcdf_read_data_1d(fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
Read in data from a variable in a netcdf file.
subroutine unpack_data_5d(fileobj, varid, varname, var_data)
subroutine unpack_data_0d(fileobj, varid, varname, var_data)
subroutine netcdf_read_data_5d(fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
Read in data from a variable in a netcdf file.
subroutine netcdf_read_data_3d(fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
Read in data from a variable in a netcdf file.
subroutine char_read_1d(fileobj, variable_name, buf, c, append_error_msg, err, varid)
Character read 1d function.
subroutine netcdf_read_data_2d(fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
Read in data from a variable in a netcdf file.
subroutine unpack_data_3d(fileobj, varid, varname, var_data)
subroutine char_read_0d(fileobj, variable_name, buf, corner, append_error_msg, err, varid)
Character read 0d function.
subroutine unpack_data_4d(fileobj, varid, varname, var_data)
subroutine netcdf_read_data_0d(fileobj, variable_name, buf, unlim_dim_level, corner, broadcast)
Read in data from a variable in a netcdf file.
subroutine netcdf_read_data_4d(fileobj, variable_name, buf, unlim_dim_level, corner, edge_lengths, broadcast)
Read in data from a variable in a netcdf file.