32 class(fmsnetcdffile_t),
intent(in) :: fileobj
33 character(len=*),
intent(in) :: variable_name
34 class(*),
intent(in) :: cdata
38 integer,
intent(in),
optional :: unlim_dim_level
40 integer,
intent(in),
optional :: corner
45 integer,
dimension(2) :: compressed_dim_index
47 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
48 if (compressed_dim_index(1) .eq. dimension_not_found)
then
49 call netcdf_write_data(fileobj, variable_name, cdata, &
50 unlim_dim_level=unlim_dim_level, corner=corner)
60 type(fmsnetcdffile_t),
intent(in) :: fileobj
61 character(len=*),
intent(in) :: variable_name
62 class(*),
intent(in) :: cdata
66 integer,
intent(in),
optional :: unlim_dim_level
68 integer,
intent(in),
optional :: corner
73 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner)
84 class(fmsnetcdffile_t),
intent(in) :: fileobj
85 character(len=*),
intent(in) :: variable_name
86 class(*),
dimension(:),
intent(in) :: cdata
90 integer,
intent(in),
optional :: unlim_dim_level
92 integer,
dimension(1),
intent(in),
optional :: corner
96 integer,
dimension(1),
intent(in),
optional :: edge_lengths
101 integer,
dimension(2) :: compressed_dim_index
104 integer,
dimension(1) :: e
106 integer(kind=i4_kind),
dimension(:),
allocatable :: buf_i4_kind
107 integer(kind=i8_kind),
dimension(:),
allocatable :: buf_i8_kind
108 real (kind=r4_kind),
dimension(:),
allocatable :: buf_r4_kind
109 real (kind=r8_kind),
dimension(:),
allocatable :: buf_r8_kind
111 character(len=200) :: append_error_msg
113 append_error_msg =
"compressed_write_1d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
115 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
117 if (compressed_dim_index(1) .eq. dimension_not_found)
then
118 call netcdf_write_data(fileobj, variable_name, cdata, &
119 unlim_dim_level=unlim_dim_level, corner=corner, &
120 edge_lengths=edge_lengths)
126 if (fileobj%is_root)
then
127 e(compressed_dim_index(1)) = sum(fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems)
129 type is (
integer(kind=i4_kind))
130 call allocate_array(buf_i4_kind, e)
131 type is (
integer(kind=i8_kind))
132 call allocate_array(buf_i8_kind, e)
133 type is (real(kind=r4_kind))
134 call allocate_array(buf_r4_kind, e)
135 type is (real(kind=r8_kind))
136 call allocate_array(buf_r8_kind, e)
138 call error(
"unsupported variable type: "//trim(append_error_msg))
142 type is (
integer(kind=i4_kind))
143 allocate(buf_i4_kind(1))
144 type is (
integer(kind=i8_kind))
145 allocate(buf_i8_kind(1))
146 type is (real(kind=r4_kind))
147 allocate(buf_r4_kind(1))
148 type is (real(kind=r8_kind))
149 allocate(buf_r8_kind(1))
151 call error(
"unsupported variable type: "//trim(append_error_msg))
157 type is (
integer(kind=i4_kind))
158 call mpp_gather(cdata,
size(cdata), buf_i4_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
160 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
161 unlim_dim_level=unlim_dim_level)
162 type is (
integer(kind=i8_kind))
163 call mpp_gather(cdata,
size(cdata), buf_i8_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
165 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
166 unlim_dim_level=unlim_dim_level)
167 type is (real(kind=r4_kind))
168 call mpp_gather(cdata,
size(cdata), buf_r4_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
170 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
171 unlim_dim_level=unlim_dim_level)
172 type is (real(kind=r8_kind))
173 call mpp_gather(cdata,
size(cdata), buf_r8_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
175 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
176 unlim_dim_level=unlim_dim_level)
178 call error(
"unsupported variable type: "//trim(append_error_msg))
180 if (
allocated(buf_i4_kind))
deallocate(buf_i4_kind)
181 if (
allocated(buf_i8_kind))
deallocate(buf_i8_kind)
182 if (
allocated(buf_r4_kind))
deallocate(buf_r4_kind)
183 if (
allocated(buf_r8_kind))
deallocate(buf_r8_kind)
192 corner, edge_lengths)
194 class(fmsnetcdffile_t),
intent(in) :: fileobj
195 character(len=*),
intent(in) :: variable_name
196 class(*),
dimension(:,:),
intent(in) :: cdata
200 integer,
intent(in),
optional :: unlim_dim_level
202 integer,
dimension(2),
intent(in),
optional :: corner
206 integer,
dimension(2),
intent(in),
optional :: edge_lengths
211 integer,
dimension(2) :: compressed_dim_index
214 integer,
dimension(2) :: c
215 integer,
dimension(2) :: e
217 integer(kind=i4_kind),
dimension(:,:),
allocatable :: buf_i4_kind
218 integer(kind=i8_kind),
dimension(:,:),
allocatable :: buf_i8_kind
219 real (kind=r4_kind),
dimension(:,:),
allocatable :: buf_r4_kind
220 real (kind=r8_kind),
dimension(:,:),
allocatable :: buf_r8_kind
222 character(len=200) :: append_error_msg
231 append_error_msg =
"compressed_write_2d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
233 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
234 if (compressed_dim_index(1) .eq. dimension_not_found)
then
235 call netcdf_write_data(fileobj, variable_name, cdata, &
236 unlim_dim_level=unlim_dim_level, corner=corner, &
237 edge_lengths=edge_lengths)
243 if (fileobj%is_root)
then
244 e(compressed_dim_index(1)) = sum(fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems)
246 type is (
integer(kind=i4_kind))
247 call allocate_array(buf_i4_kind, e)
248 type is (
integer(kind=i8_kind))
249 call allocate_array(buf_i8_kind, e)
250 type is (real(kind=r4_kind))
251 call allocate_array(buf_r4_kind, e)
252 type is (real(kind=r8_kind))
253 call allocate_array(buf_r8_kind, e)
255 call error(
"unsupported variable type: "//trim(append_error_msg))
260 index = findloc(fileobj%pelist,
mpp_pe())
262 c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(index(1))
263 e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(index(1))
265 if (compressed_dim_index(1) .eq. 1)
then
266 is = c(compressed_dim_index(1))
267 ie = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
271 js = c(compressed_dim_index(1))
272 je = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
278 type is (
integer(kind=i4_kind))
279 call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_i4_kind, fileobj%is_root)
280 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
281 unlim_dim_level=unlim_dim_level)
282 type is (
integer(kind=i8_kind))
283 call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_i8_kind, fileobj%is_root)
284 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
285 unlim_dim_level=unlim_dim_level)
286 type is (real(kind=r4_kind))
287 call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_r4_kind, fileobj%is_root)
288 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
289 unlim_dim_level=unlim_dim_level)
290 type is (real(kind=r8_kind))
291 call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_r8_kind, fileobj%is_root)
292 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
293 unlim_dim_level=unlim_dim_level)
295 call error(
"unsupported variable type: "//trim(append_error_msg))
298 if (fileobj%is_root)
then
299 if (
allocated(buf_i4_kind))
deallocate(buf_i4_kind)
300 if (
allocated(buf_i8_kind))
deallocate(buf_i8_kind)
301 if (
allocated(buf_r4_kind))
deallocate(buf_r4_kind)
302 if (
allocated(buf_r8_kind))
deallocate(buf_r8_kind)
312 corner, edge_lengths)
314 class(fmsnetcdffile_t),
intent(in) :: fileobj
315 character(len=*),
intent(in) :: variable_name
316 class(*),
contiguous,
intent(in),
target :: cdata(:,:,:)
320 integer,
intent(in),
optional :: unlim_dim_level
322 integer,
dimension(3),
intent(in),
optional :: corner
326 integer,
dimension(3),
intent(in),
optional :: edge_lengths
331 integer,
dimension(2) :: compressed_dim_index
334 integer,
dimension(3) :: c
335 integer,
dimension(3) :: e
337 integer(kind=i4_kind),
dimension(:,:,:),
allocatable :: buf_i4_kind
338 integer(kind=i8_kind),
dimension(:,:,:),
allocatable :: buf_i8_kind
339 real (kind=r4_kind),
dimension(:,:,:),
allocatable :: buf_r4_kind
340 real (kind=r8_kind),
dimension(:,:,:),
allocatable :: buf_r8_kind
342 character(len=200) :: append_error_msg
343 class(*),
pointer :: cdata_dummy(:,:,:,:)
352 append_error_msg =
"compressed_write_3d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
354 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
355 if (compressed_dim_index(1) .eq. dimension_not_found)
then
356 call netcdf_write_data(fileobj, variable_name, cdata, &
357 unlim_dim_level=unlim_dim_level, corner=corner, &
358 edge_lengths=edge_lengths)
360 else if (compressed_dim_index(1) .eq. 3)
then
361 cdata_dummy(1:
size(cdata,1), 1:
size(cdata,2), 1:
size(cdata,3), 1:1) => cdata(:,:,:)
367 if (fileobj%is_root)
then
368 e(compressed_dim_index(1)) = sum(fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems)
370 type is (
integer(kind=i4_kind))
371 call allocate_array(buf_i4_kind, e)
372 type is (
integer(kind=i8_kind))
373 call allocate_array(buf_i8_kind, e)
374 type is (real(kind=r4_kind))
375 call allocate_array(buf_r4_kind, e)
376 type is (real(kind=r8_kind))
377 call allocate_array(buf_r8_kind, e)
379 call error(
"unsupported variable type: "//trim(append_error_msg))
384 index = findloc(fileobj%pelist,
mpp_pe())
386 c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(index(1))
387 e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(index(1))
389 if (compressed_dim_index(1) .eq. 1)
then
390 is = c(compressed_dim_index(1))
391 ie = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
395 js = c(compressed_dim_index(1))
396 je = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
402 type is (
integer(kind=i4_kind))
403 call mpp_gather(is, ie, js, je,
size(cdata,3), &
404 fileobj%pelist, cdata, buf_i4_kind, fileobj%is_root)
405 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
406 unlim_dim_level=unlim_dim_level)
407 type is (
integer(kind=i8_kind))
408 call mpp_gather(is, ie, js, je,
size(cdata,3), &
409 fileobj%pelist, cdata, buf_i8_kind, fileobj%is_root)
410 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
411 unlim_dim_level=unlim_dim_level)
412 type is (real(kind=r4_kind))
413 call mpp_gather(is, ie, js, je,
size(cdata,3), &
414 fileobj%pelist, cdata, buf_r4_kind, fileobj%is_root)
415 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
416 unlim_dim_level=unlim_dim_level)
417 type is (real(kind=r8_kind))
418 call mpp_gather(is, ie, js, je,
size(cdata,3), &
419 fileobj%pelist, cdata, buf_r8_kind, fileobj%is_root)
420 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
421 unlim_dim_level=unlim_dim_level)
423 call error(
"unsupported variable type: "//trim(append_error_msg))
426 if (fileobj%is_root)
then
427 if (
allocated(buf_i4_kind))
deallocate(buf_i4_kind)
428 if (
allocated(buf_i8_kind))
deallocate(buf_i8_kind)
429 if (
allocated(buf_r4_kind))
deallocate(buf_r4_kind)
430 if (
allocated(buf_r8_kind))
deallocate(buf_r8_kind)
440 corner, edge_lengths)
442 class(fmsnetcdffile_t),
intent(in) :: fileobj
443 character(len=*),
intent(in) :: variable_name
444 class(*),
dimension(:,:,:,:),
intent(in) :: cdata
448 integer,
intent(in),
optional :: unlim_dim_level
450 integer,
dimension(4),
intent(in),
optional :: corner
454 integer,
dimension(4),
intent(in),
optional :: edge_lengths
459 integer,
dimension(2) :: compressed_dim_index
460 integer,
dimension(4) :: c
461 integer,
dimension(4) :: e
463 integer(kind=i4_kind),
dimension(:,:,:,:),
allocatable :: buf_i4_kind
464 integer(kind=i8_kind),
dimension(:,:,:,:),
allocatable :: buf_i8_kind
465 real(kind=r4_kind),
dimension(:,:,:,:),
allocatable :: buf_r4_kind
466 real(kind=r8_kind),
dimension(:,:,:,:),
allocatable :: buf_r8_kind
467 character(len=200) :: append_error_msg
469 append_error_msg =
"compressed_write_4d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
471 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
472 if (compressed_dim_index(1) .eq. dimension_not_found)
then
473 call netcdf_write_data(fileobj, variable_name, cdata, &
474 unlim_dim_level=unlim_dim_level, corner=corner, &
475 edge_lengths=edge_lengths)
480 if (fileobj%is_root)
then
483 do i = 1,
size(fileobj%pelist)
484 c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(i)
485 e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(i)
487 call netcdf_write_data(fileobj, variable_name, cdata, &
488 unlim_dim_level=unlim_dim_level, corner=c, &
492 type is (
integer(kind=i4_kind))
493 call allocate_array(buf_i4_kind, e)
494 call mpp_recv(buf_i4_kind,
size(buf_i4_kind), fileobj%pelist(i), block=.true.)
495 call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
496 unlim_dim_level=unlim_dim_level, corner=c, &
498 deallocate(buf_i4_kind)
499 type is (
integer(kind=i8_kind))
500 call allocate_array(buf_i8_kind, e)
501 call mpp_recv(buf_i8_kind,
size(buf_i8_kind), fileobj%pelist(i), block=.true.)
502 call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
503 unlim_dim_level=unlim_dim_level, corner=c, &
505 deallocate(buf_i8_kind)
506 type is (real(kind=r4_kind))
507 call allocate_array(buf_r4_kind, e)
508 call mpp_recv(buf_r4_kind,
size(buf_r4_kind), fileobj%pelist(i), block=.true.)
509 call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
510 unlim_dim_level=unlim_dim_level, corner=c, &
512 deallocate(buf_r4_kind)
513 type is (real(kind=r8_kind))
514 call allocate_array(buf_r8_kind, e)
515 call mpp_recv(buf_r8_kind,
size(buf_r8_kind), fileobj%pelist(i), block=.true.)
516 call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
517 unlim_dim_level=unlim_dim_level, corner=c, &
519 deallocate(buf_r8_kind)
521 call error(
"unsupported variable type: "//trim(append_error_msg))
527 type is (
integer(kind=i4_kind))
528 call mpp_send(cdata,
size(cdata), fileobj%io_root)
529 type is (
integer(kind=i8_kind))
530 call mpp_send(cdata,
size(cdata), fileobj%io_root)
531 type is (real(kind=r4_kind))
532 call mpp_send(cdata,
size(cdata), fileobj%io_root)
533 type is (real(kind=r8_kind))
534 call mpp_send(cdata,
size(cdata), fileobj%io_root)
536 call error(
"unsupported variable type: "//trim(append_error_msg))
548 corner, edge_lengths)
550 class(fmsnetcdffile_t),
intent(in) :: fileobj
551 character(len=*),
intent(in) :: variable_name
552 class(*),
dimension(:,:,:,:,:),
intent(in) :: cdata
556 integer,
intent(in),
optional :: unlim_dim_level
558 integer,
dimension(5),
intent(in),
optional :: corner
562 integer,
dimension(5),
intent(in),
optional :: edge_lengths
567 integer,
dimension(2) :: compressed_dim_index
568 integer,
dimension(5) :: c
569 integer,
dimension(5) :: e
571 integer(kind=i4_kind),
dimension(:,:,:,:,:),
allocatable :: buf_i4_kind
572 integer(kind=i8_kind),
dimension(:,:,:,:,:),
allocatable :: buf_i8_kind
573 real(kind=r4_kind),
dimension(:,:,:,:,:),
allocatable :: buf_r4_kind
574 real(kind=r8_kind),
dimension(:,:,:,:,:),
allocatable :: buf_r8_kind
575 character(len=200) :: append_error_msg
577 append_error_msg =
"compressed_write_5d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
579 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
580 if (compressed_dim_index(1) .eq. dimension_not_found)
then
581 call netcdf_write_data(fileobj, variable_name, cdata, &
582 unlim_dim_level=unlim_dim_level, corner=corner, &
583 edge_lengths=edge_lengths)
588 if (fileobj%is_root)
then
591 do i = 1,
size(fileobj%pelist)
592 c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(i)
593 e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(i)
595 call netcdf_write_data(fileobj, variable_name, cdata, &
596 unlim_dim_level=unlim_dim_level, corner=c, &
600 type is (
integer(kind=i4_kind))
601 call allocate_array(buf_i4_kind, e)
602 call mpp_recv(buf_i4_kind,
size(buf_i4_kind), fileobj%pelist(i), block=.true.)
603 call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
604 unlim_dim_level=unlim_dim_level, corner=c, &
606 deallocate(buf_i4_kind)
607 type is (
integer(kind=i8_kind))
608 call allocate_array(buf_i8_kind, e)
609 call mpp_recv(buf_i8_kind,
size(buf_i8_kind), fileobj%pelist(i), block=.true.)
610 call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
611 unlim_dim_level=unlim_dim_level, corner=c, &
613 deallocate(buf_i8_kind)
614 type is (real(kind=r4_kind))
615 call allocate_array(buf_r4_kind, e)
616 call mpp_recv(buf_r4_kind,
size(buf_r4_kind), fileobj%pelist(i), block=.true.)
617 call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
618 unlim_dim_level=unlim_dim_level, corner=c, &
620 deallocate(buf_r4_kind)
621 type is (real(kind=r8_kind))
622 call allocate_array(buf_r8_kind, e)
623 call mpp_recv(buf_r8_kind,
size(buf_r8_kind), fileobj%pelist(i), block=.true.)
624 call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
625 unlim_dim_level=unlim_dim_level, corner=c, &
627 deallocate(buf_r8_kind)
629 call error(
"unsupported variable type: "//trim(append_error_msg))
635 type is (
integer(kind=i4_kind))
636 call mpp_send(cdata,
size(cdata), fileobj%io_root)
637 type is (
integer(kind=i8_kind))
638 call mpp_send(cdata,
size(cdata), fileobj%io_root)
639 type is (real(kind=r4_kind))
640 call mpp_send(cdata,
size(cdata), fileobj%io_root)
641 type is (real(kind=r8_kind))
642 call mpp_send(cdata,
size(cdata), fileobj%io_root)
644 call error(
"unsupported variable type: "//trim(append_error_msg))
653 corner, edge_lengths)
655 type(fmsnetcdffile_t),
intent(in) :: fileobj
656 character(len=*),
intent(in) :: variable_name
657 class(*),
dimension(:),
intent(in) :: cdata
661 integer,
intent(in),
optional :: unlim_dim_level
663 integer,
dimension(1),
intent(in),
optional :: corner
667 integer,
dimension(1),
intent(in),
optional :: edge_lengths
672 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
673 edge_lengths=edge_lengths)
680 corner, edge_lengths)
682 type(fmsnetcdffile_t),
intent(in) :: fileobj
683 character(len=*),
intent(in) :: variable_name
684 class(*),
dimension(:,:),
intent(in) :: cdata
688 integer,
intent(in),
optional :: unlim_dim_level
690 integer,
dimension(2),
intent(in),
optional :: corner
694 integer,
dimension(2),
intent(in),
optional :: edge_lengths
699 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
700 edge_lengths=edge_lengths)
706 corner, edge_lengths)
708 type(fmsnetcdffile_t),
intent(in) :: fileobj
709 character(len=*),
intent(in) :: variable_name
710 class(*),
dimension(:,:,:),
intent(in) :: cdata
714 integer,
intent(in),
optional :: unlim_dim_level
716 integer,
dimension(3),
intent(in),
optional :: corner
720 integer,
dimension(3),
intent(in),
optional :: edge_lengths
725 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
726 edge_lengths=edge_lengths)
732 corner, edge_lengths)
734 type(fmsnetcdffile_t),
intent(in) :: fileobj
735 character(len=*),
intent(in) :: variable_name
736 class(*),
dimension(:,:,:,:),
intent(in) :: cdata
740 integer,
intent(in),
optional :: unlim_dim_level
742 integer,
dimension(4),
intent(in),
optional :: corner
746 integer,
dimension(4),
intent(in),
optional :: edge_lengths
751 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
752 edge_lengths=edge_lengths)
758 corner, edge_lengths)
760 type(fmsnetcdffile_t),
intent(in) :: fileobj
761 character(len=*),
intent(in) :: variable_name
762 class(*),
dimension(:,:,:,:,:),
intent(in) :: cdata
766 integer,
intent(in),
optional :: unlim_dim_level
768 integer,
dimension(5),
intent(in),
optional :: corner
772 integer,
dimension(5),
intent(in),
optional :: edge_lengths
777 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
778 edge_lengths=edge_lengths)
subroutine compressed_write_1d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine compressed_write_3d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_5d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine compressed_write_0d(fileobj, variable_name, cdata, unlim_dim_level, corner)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_5d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_0d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner)
Wrapper to distinguish interfaces.
subroutine compressed_write_4d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_2d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine compressed_write_1d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_4d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine compressed_write_2d(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
For variables without a compressed dimension, this routine simply wraps netcdf_write data....
subroutine compressed_write_3d_wrap(fileobj, variable_name, cdata, unlim_dim_level, corner, edge_lengths)
Wrapper to distinguish interfaces.
subroutine mpp_sync_self(pelist, check, request, msg_size, msg_type)
This is to check if current PE's outstanding puts are complete but we can't use shmem_fence because w...
integer function mpp_pe()
Returns processor ID.