27 subroutine char_read_0d(fileobj, variable_name, buf, corner, append_error_msg, err, varid)
28 class(fmsnetcdffile_t),
intent(in) :: fileobj
29 character(len=*),
intent(in) :: variable_name
30 character(len=*),
intent(inout) :: buf
31 integer,
intent(in),
optional :: corner
33 character(len=200),
intent(in):: append_error_msg
34 integer,
intent(inout) :: err
35 integer,
intent(in) :: varid
37 integer,
dimension(2) :: start
39 character,
dimension(:),
allocatable :: charbuf
40 integer,
dimension(:),
allocatable :: dimsizes
44 ndims = get_variable_num_dimensions(fileobj, variable_name, broadcast=.false.)
45 allocate(dimsizes(ndims))
46 call get_variable_size(fileobj, variable_name, dimsizes, broadcast=.false.)
47 allocate(charbuf(dimsizes(1)))
49 if (ndims .eq. 2)
then
50 if (
present(corner))
then
54 elseif (ndims .gt. 2)
then
55 call error(
"Only scalar and 1d string values are currently supported: "//trim(append_error_msg))
57 err = nf90_get_var(fileobj%ncid, varid, charbuf, start=start, count=dimsizes)
58 if (len(buf) .lt. dimsizes(1))
then
59 call error(
"character buffer is too small; increase length: "//trim(append_error_msg))
63 if (charbuf(i) .eq. char(0))
then
73 subroutine char_read_1d(fileobj, variable_name, buf, c, append_error_msg, err, varid)
75 class(fmsnetcdffile_t),
intent(in) :: fileobj
76 character(len=*),
intent(in) :: variable_name
77 character(len=*),
dimension(:),
intent(inout) :: buf
79 integer,
dimension(2),
intent(in) :: c
80 character(len=200),
intent(in) :: append_error_msg
81 integer,
intent(inout) :: err
82 integer,
intent(in) :: varid
85 integer,
dimension(2) :: start
86 integer,
dimension(2) :: dimsizes
87 character,
dimension(:,:),
allocatable :: charbuf
88 character(len=1024) :: sbuf
92 ndims = get_variable_num_dimensions(fileobj, variable_name, broadcast=.false.)
93 if (ndims .ne. 2)
then
94 call error(trim(variable_name)//
"must be 2d dimensional in netcdf file.")
98 call get_variable_size(fileobj, variable_name, dimsizes, .false.)
99 dimsizes(2) = dimsizes(2) - start(2) + 1
100 call allocate_array(charbuf, dimsizes, initialize=.true.)
101 err = nf90_get_var(fileobj%ncid, varid, charbuf, start=start, count=dimsizes)
102 if (len(buf(1)) .lt. dimsizes(1))
then
103 call error(
"character buffer is too small; increase length: "//trim(append_error_msg))
105 if (
size(buf) .lt. dimsizes(2))
then
106 call error(
"incorrect buffer array size:: "//trim(append_error_msg))
108 do i = start(2), start(2)+dimsizes(2)-1
110 do j = 1, dimsizes(1)
111 if (charbuf(j, i-start(2)+1) .eq. char(0))
then
114 sbuf(j:j) = charbuf(j, i-start(2)+1)
116 call string_copy(buf(i-start(2)+1), sbuf)
126 class(fmsnetcdffile_t),
intent(in) :: fileobj
127 character(len=*),
intent(in) :: variable_name
128 class(*),
intent(inout) :: buf
129 integer,
intent(in),
optional :: unlim_dim_level
130 integer,
intent(in),
optional :: corner
132 logical,
intent(in),
optional :: broadcast
142 integer :: unlim_dim_index
143 integer,
dimension(1) :: c
144 character(len=1024),
dimension(1) :: buf1d
145 character(len=200) :: append_error_msg
147 append_error_msg =
"netcdf_read_data_0d: file:"//trim(fileobj%path)//
"- variable:"//trim(variable_name)
149 if (
present(broadcast))
then
155 if (
present(unlim_dim_level))
then
156 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
158 if (unlim_dim_index .ne. 1)
then
159 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
161 c(unlim_dim_index) = unlim_dim_level
163 if (fileobj%is_root)
then
164 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
166 type is (
integer(kind=i4_kind))
167 err = nf90_get_var(fileobj%ncid, varid, buf, start=c)
168 type is (
integer(kind=i8_kind))
169 err = nf90_get_var(fileobj%ncid, varid, buf, start=c)
170 type is (real(kind=r4_kind))
171 err = nf90_get_var(fileobj%ncid, varid, buf, start=c)
172 type is (real(kind=r8_kind))
173 err = nf90_get_var(fileobj%ncid, varid, buf, start=c)
174 type is (
character(len=*))
175 call char_read_0d(fileobj, variable_name, buf, corner, append_error_msg, err, varid)
177 call error(
"Unsupported variable type: "//trim(append_error_msg))
179 call check_netcdf_code(err, append_error_msg)
184 type is (
integer(kind=i4_kind))
185 call mpp_broadcast(buf, fileobj%io_root, pelist=fileobj%pelist)
186 type is (
integer(kind=i8_kind))
187 call mpp_broadcast(buf, fileobj%io_root, pelist=fileobj%pelist)
188 type is (real(kind=r4_kind))
189 call mpp_broadcast(buf, fileobj%io_root, pelist=fileobj%pelist)
190 type is (real(kind=r8_kind))
191 call mpp_broadcast(buf, fileobj%io_root, pelist=fileobj%pelist)
192 type is (
character(len=*))
193 call string_copy(buf1d(1), buf)
194 call mpp_broadcast(buf1d, len(buf1d(1)), fileobj%io_root, pelist=fileobj%pelist)
195 call string_copy(buf, buf1d(1))
197 call error(
"Unsupported variable type: "//trim(append_error_msg))
205 corner, edge_lengths, broadcast)
207 class(fmsnetcdffile_t),
intent(in) :: fileobj
208 character(len=*),
intent(in) :: variable_name
209 class(*),
dimension(:),
intent(inout) :: buf
211 integer,
intent(in),
optional :: unlim_dim_level
213 integer,
dimension(:),
intent(in),
optional :: corner
217 integer,
dimension(:),
intent(in),
optional :: edge_lengths
221 logical,
intent(in),
optional :: broadcast
231 integer :: unlim_dim_index
232 integer,
dimension(2) :: c
233 integer,
dimension(2) :: e
234 character(len=200) :: append_error_msg
236 append_error_msg =
"netcdf_read_data_1d: file:"//trim(fileobj%path)//
"- variable:"//trim(variable_name)
238 if (
present(broadcast))
then
244 if (
present(corner))
then
248 if (
present(edge_lengths))
then
249 e(1) = edge_lengths(1)
253 if (
present(unlim_dim_level))
then
254 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
256 if (unlim_dim_index .ne. 2)
then
257 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
259 c(unlim_dim_index) = unlim_dim_level
261 if (fileobj%is_root)
then
262 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
264 type is (
integer(kind=i4_kind))
265 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
266 type is (
integer(kind=i8_kind))
267 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
268 type is (real(kind=r4_kind))
269 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
270 type is (real(kind=r8_kind))
271 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
272 type is (
character(len=*))
273 call char_read_1d(fileobj, variable_name, buf, c, append_error_msg, err, varid)
275 call error(
"Unsupported variable type: "//trim(append_error_msg))
277 call check_netcdf_code(err, append_error_msg)
282 type is (
integer(kind=i4_kind))
283 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
284 type is (
integer(kind=i8_kind))
285 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
286 type is (real(kind=r4_kind))
287 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
288 type is (real(kind=r8_kind))
289 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
290 type is (
character(len=*))
291 call mpp_broadcast(buf, len(buf(1)), fileobj%io_root, pelist=fileobj%pelist)
293 call error(
"Unsupported variable type: "//trim(append_error_msg))
301 corner, edge_lengths, broadcast)
303 class(fmsnetcdffile_t),
intent(in) :: fileobj
304 character(len=*),
intent(in) :: variable_name
305 class(*),
dimension(:,:),
intent(inout) :: buf
307 integer,
intent(in),
optional :: unlim_dim_level
309 integer,
dimension(:),
intent(in),
optional :: corner
313 integer,
dimension(:),
intent(in),
optional :: edge_lengths
317 logical,
intent(in),
optional :: broadcast
327 integer :: unlim_dim_index
328 integer,
dimension(3) :: c
329 integer,
dimension(3) :: e
330 character(len=200) :: append_error_msg
332 append_error_msg =
"netcdf_read_data_2d: file:"//trim(fileobj%path)//
"- variable:"//trim(variable_name)
334 if (
present(broadcast))
then
340 if (
present(corner))
then
344 if (
present(edge_lengths))
then
345 e(1:2) = edge_lengths(1:2)
349 if (
present(unlim_dim_level))
then
350 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
352 if (unlim_dim_index .ne. 3)
then
353 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
355 c(unlim_dim_index) = unlim_dim_level
357 if(fileobj%use_collective)
then
358 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
361 err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
362 call check_netcdf_code(err, append_error_msg)
364 type is (
integer(kind=i4_kind))
365 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
366 type is (
integer(kind=i8_kind))
367 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
368 type is (real(kind=r4_kind))
369 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
370 type is (real(kind=r8_kind))
371 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
373 call error(
"Unsupported variable type: "//trim(append_error_msg))
375 call check_netcdf_code(err, append_error_msg)
378 if (fileobj%is_root)
then
379 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
381 type is (
integer(kind=i4_kind))
382 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
383 type is (
integer(kind=i8_kind))
384 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
385 type is (real(kind=r4_kind))
386 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
387 type is (real(kind=r8_kind))
388 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
390 call error(
"Unsupported variable type: "//trim(append_error_msg))
392 call check_netcdf_code(err, append_error_msg)
397 type is (
integer(kind=i4_kind))
398 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
399 type is (
integer(kind=i8_kind))
400 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
401 type is (real(kind=r4_kind))
402 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
403 type is (real(kind=r8_kind))
404 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
406 call error(
"Unsupported variable type: "//trim(append_error_msg))
415 corner, edge_lengths, broadcast)
417 class(fmsnetcdffile_t),
intent(in) :: fileobj
418 character(len=*),
intent(in) :: variable_name
419 class(*),
dimension(:,:,:),
intent(inout) :: buf
421 integer,
intent(in),
optional :: unlim_dim_level
423 integer,
dimension(:),
intent(in),
optional :: corner
427 integer,
dimension(:),
intent(in),
optional :: edge_lengths
431 logical,
intent(in),
optional :: broadcast
441 integer :: unlim_dim_index
442 integer,
dimension(4) :: c
443 integer,
dimension(4) :: e
444 character(len=200) :: append_error_msg
446 append_error_msg =
"netcdf_read_data_3d: file:"//trim(fileobj%path)//
"- variable:"//trim(variable_name)
448 if (
present(broadcast))
then
454 if (
present(corner))
then
458 if (
present(edge_lengths))
then
459 e(1:3) = edge_lengths(1:3)
463 if (
present(unlim_dim_level))
then
464 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
466 if (unlim_dim_index .ne. 4)
then
467 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
469 c(unlim_dim_index) = unlim_dim_level
471 if(fileobj%use_collective)
then
472 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
475 err = nf90_var_par_access(fileobj%ncid, varid, nf90_collective)
476 call check_netcdf_code(err, append_error_msg)
478 type is (
integer(kind=i4_kind))
479 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
480 type is (
integer(kind=i8_kind))
481 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
482 type is (real(kind=r4_kind))
483 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
484 type is (real(kind=r8_kind))
485 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
487 call error(
"Unsupported variable type: "//trim(append_error_msg))
489 call check_netcdf_code(err, append_error_msg)
492 if (fileobj%is_root)
then
493 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
495 type is (
integer(kind=i4_kind))
496 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
497 type is (
integer(kind=i8_kind))
498 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
499 type is (real(kind=r4_kind))
500 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
501 type is (real(kind=r8_kind))
502 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
504 call error(
"Unsupported variable type: "//trim(append_error_msg))
506 call check_netcdf_code(err, append_error_msg)
511 type is (
integer(kind=i4_kind))
512 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
513 type is (
integer(kind=i8_kind))
514 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
515 type is (real(kind=r4_kind))
516 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
517 type is (real(kind=r8_kind))
518 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
520 call error(
"Unsupported variable type: "//trim(append_error_msg))
529 corner, edge_lengths, broadcast)
531 class(fmsnetcdffile_t),
intent(in) :: fileobj
532 character(len=*),
intent(in) :: variable_name
533 class(*),
dimension(:,:,:,:),
intent(inout) :: buf
535 integer,
intent(in),
optional :: unlim_dim_level
537 integer,
dimension(:),
intent(in),
optional :: corner
541 integer,
dimension(:),
intent(in),
optional :: edge_lengths
545 logical,
intent(in),
optional :: broadcast
555 integer :: unlim_dim_index
556 integer,
dimension(5) :: c
557 integer,
dimension(5) :: e
558 character(len=200) :: append_error_msg
560 append_error_msg =
"netcdf_read_data_4d: file:"//trim(fileobj%path)//
"- variable:"//trim(variable_name)
562 if (
present(broadcast))
then
568 if (
present(corner))
then
572 if (
present(edge_lengths))
then
573 e(1:4) = edge_lengths(1:4)
577 if (
present(unlim_dim_level))
then
578 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
580 if (unlim_dim_index .ne. 5)
then
581 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
583 c(unlim_dim_index) = unlim_dim_level
585 if (fileobj%is_root)
then
586 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
588 type is (
integer(kind=i4_kind))
589 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
590 type is (
integer(kind=i8_kind))
591 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
592 type is (real(kind=r4_kind))
593 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
594 type is (real(kind=r8_kind))
595 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
597 call error(
"Unsupported variable type: "//trim(append_error_msg))
599 call check_netcdf_code(err, append_error_msg)
604 type is (
integer(kind=i4_kind))
605 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
606 type is (
integer(kind=i8_kind))
607 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
608 type is (real(kind=r4_kind))
609 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
610 type is (real(kind=r8_kind))
611 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
613 call error(
"Unsupported variable type: "//trim(append_error_msg))
621 corner, edge_lengths, broadcast)
623 class(fmsnetcdffile_t),
intent(in) :: fileobj
624 character(len=*),
intent(in) :: variable_name
625 class(*),
dimension(:,:,:,:,:),
intent(inout) :: buf
627 integer,
intent(in),
optional :: unlim_dim_level
629 integer,
dimension(:),
intent(in),
optional :: corner
633 integer,
dimension(:),
intent(in),
optional :: edge_lengths
637 logical,
intent(in),
optional :: broadcast
647 integer :: unlim_dim_index
648 integer,
dimension(6) :: c
649 integer,
dimension(6) :: e
650 character(len=200) :: append_error_msg
652 append_error_msg =
"netcdf_read_data_5d: file:"//trim(fileobj%path)//
"- variable:"//trim(variable_name)
654 if (
present(broadcast))
then
660 if (
present(corner))
then
664 if (
present(edge_lengths))
then
665 e(1:5) = edge_lengths(1:5)
669 if (
present(unlim_dim_level))
then
670 unlim_dim_index = get_variable_unlimited_dimension_index(fileobj, variable_name, &
672 if (unlim_dim_index .ne. 6)
then
673 call error(
"unlimited dimension must be the slowest varying dimension: "//trim(append_error_msg))
675 c(unlim_dim_index) = unlim_dim_level
677 if (fileobj%is_root)
then
678 varid = get_variable_id(fileobj%ncid, trim(variable_name), msg=append_error_msg)
680 type is (
integer(kind=i4_kind))
681 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
682 type is (
integer(kind=i8_kind))
683 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
684 type is (real(kind=r4_kind))
685 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
686 type is (real(kind=r8_kind))
687 err = nf90_get_var(fileobj%ncid, varid, buf, start=c, count=e)
689 call error(
"Unsupported variable type: "//trim(append_error_msg))
691 call check_netcdf_code(err, append_error_msg)
696 type is (
integer(kind=i4_kind))
697 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
698 type is (
integer(kind=i8_kind))
699 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
700 type is (real(kind=r4_kind))
701 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
702 type is (real(kind=r8_kind))
703 call mpp_broadcast(buf,
size(buf), fileobj%io_root, pelist=fileobj%pelist)
705 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.