31 class(fmsnetcdffile_t),
intent(in) :: fileobj
32 character(len=*),
intent(in) :: variable_name
33 class(*),
intent(in) :: cdata
37 integer,
intent(in),
optional :: unlim_dim_level
39 integer,
intent(in),
optional :: corner
44 integer,
dimension(2) :: compressed_dim_index
46 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
47 if (compressed_dim_index(1) .eq. dimension_not_found)
then
48 call netcdf_write_data(fileobj, variable_name, cdata, &
49 unlim_dim_level=unlim_dim_level, corner=corner)
59 type(fmsnetcdffile_t),
intent(in) :: fileobj
60 character(len=*),
intent(in) :: variable_name
61 class(*),
intent(in) :: cdata
65 integer,
intent(in),
optional :: unlim_dim_level
67 integer,
intent(in),
optional :: corner
72 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner)
83 class(fmsnetcdffile_t),
intent(in) :: fileobj
84 character(len=*),
intent(in) :: variable_name
85 class(*),
dimension(:),
intent(in) :: cdata
89 integer,
intent(in),
optional :: unlim_dim_level
91 integer,
dimension(1),
intent(in),
optional :: corner
95 integer,
dimension(1),
intent(in),
optional :: edge_lengths
100 integer,
dimension(2) :: compressed_dim_index
103 integer,
dimension(1) :: e
105 integer(kind=i4_kind),
dimension(:),
allocatable :: buf_i4_kind
106 integer(kind=i8_kind),
dimension(:),
allocatable :: buf_i8_kind
107 real (kind=r4_kind),
dimension(:),
allocatable :: buf_r4_kind
108 real (kind=r8_kind),
dimension(:),
allocatable :: buf_r8_kind
110 character(len=200) :: append_error_msg
112 append_error_msg =
"compressed_write_1d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
114 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
116 if (compressed_dim_index(1) .eq. dimension_not_found)
then
117 call netcdf_write_data(fileobj, variable_name, cdata, &
118 unlim_dim_level=unlim_dim_level, corner=corner, &
119 edge_lengths=edge_lengths)
125 if (fileobj%is_root)
then
126 e(compressed_dim_index(1)) = sum(fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems)
128 type is (
integer(kind=i4_kind))
129 call allocate_array(buf_i4_kind, e)
130 type is (
integer(kind=i8_kind))
131 call allocate_array(buf_i8_kind, e)
132 type is (real(kind=r4_kind))
133 call allocate_array(buf_r4_kind, e)
134 type is (real(kind=r8_kind))
135 call allocate_array(buf_r8_kind, e)
137 call error(
"unsupported variable type: "//trim(append_error_msg))
141 type is (
integer(kind=i4_kind))
142 allocate(buf_i4_kind(1))
143 type is (
integer(kind=i8_kind))
144 allocate(buf_i8_kind(1))
145 type is (real(kind=r4_kind))
146 allocate(buf_r4_kind(1))
147 type is (real(kind=r8_kind))
148 allocate(buf_r8_kind(1))
150 call error(
"unsupported variable type: "//trim(append_error_msg))
156 type is (
integer(kind=i4_kind))
157 call mpp_gather(cdata,
size(cdata), buf_i4_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
159 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
160 unlim_dim_level=unlim_dim_level)
161 type is (
integer(kind=i8_kind))
162 call mpp_gather(cdata,
size(cdata), buf_i8_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
164 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
165 unlim_dim_level=unlim_dim_level)
166 type is (real(kind=r4_kind))
167 call mpp_gather(cdata,
size(cdata), buf_r4_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
169 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
170 unlim_dim_level=unlim_dim_level)
171 type is (real(kind=r8_kind))
172 call mpp_gather(cdata,
size(cdata), buf_r8_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
174 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
175 unlim_dim_level=unlim_dim_level)
177 call error(
"unsupported variable type: "//trim(append_error_msg))
179 if (
allocated(buf_i4_kind))
deallocate(buf_i4_kind)
180 if (
allocated(buf_i8_kind))
deallocate(buf_i8_kind)
181 if (
allocated(buf_r4_kind))
deallocate(buf_r4_kind)
182 if (
allocated(buf_r8_kind))
deallocate(buf_r8_kind)
191 corner, edge_lengths)
193 class(fmsnetcdffile_t),
intent(in) :: fileobj
194 character(len=*),
intent(in) :: variable_name
195 class(*),
dimension(:,:),
intent(in) :: cdata
199 integer,
intent(in),
optional :: unlim_dim_level
201 integer,
dimension(2),
intent(in),
optional :: corner
205 integer,
dimension(2),
intent(in),
optional :: edge_lengths
210 integer,
dimension(2) :: compressed_dim_index
213 integer,
dimension(2) :: c
214 integer,
dimension(2) :: e
216 integer(kind=i4_kind),
dimension(:,:),
allocatable :: buf_i4_kind
217 integer(kind=i8_kind),
dimension(:,:),
allocatable :: buf_i8_kind
218 real (kind=r4_kind),
dimension(:,:),
allocatable :: buf_r4_kind
219 real (kind=r8_kind),
dimension(:,:),
allocatable :: buf_r8_kind
221 character(len=200) :: append_error_msg
230 append_error_msg =
"compressed_write_2d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
232 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
233 if (compressed_dim_index(1) .eq. dimension_not_found)
then
234 call netcdf_write_data(fileobj, variable_name, cdata, &
235 unlim_dim_level=unlim_dim_level, corner=corner, &
236 edge_lengths=edge_lengths)
242 if (fileobj%is_root)
then
243 e(compressed_dim_index(1)) = sum(fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems)
245 type is (
integer(kind=i4_kind))
246 call allocate_array(buf_i4_kind, e)
247 type is (
integer(kind=i8_kind))
248 call allocate_array(buf_i8_kind, e)
249 type is (real(kind=r4_kind))
250 call allocate_array(buf_r4_kind, e)
251 type is (real(kind=r8_kind))
252 call allocate_array(buf_r8_kind, e)
254 call error(
"unsupported variable type: "//trim(append_error_msg))
259 index = findloc(fileobj%pelist,
mpp_pe())
261 c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(index(1))
262 e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(index(1))
264 if (compressed_dim_index(1) .eq. 1)
then
265 is = c(compressed_dim_index(1))
266 ie = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
270 js = c(compressed_dim_index(1))
271 je = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
277 type is (
integer(kind=i4_kind))
278 call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_i4_kind, fileobj%is_root)
279 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
280 unlim_dim_level=unlim_dim_level)
281 type is (
integer(kind=i8_kind))
282 call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_i8_kind, fileobj%is_root)
283 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
284 unlim_dim_level=unlim_dim_level)
285 type is (real(kind=r4_kind))
286 call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_r4_kind, fileobj%is_root)
287 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
288 unlim_dim_level=unlim_dim_level)
289 type is (real(kind=r8_kind))
290 call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_r8_kind, fileobj%is_root)
291 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
292 unlim_dim_level=unlim_dim_level)
294 call error(
"unsupported variable type: "//trim(append_error_msg))
297 if (fileobj%is_root)
then
298 if (
allocated(buf_i4_kind))
deallocate(buf_i4_kind)
299 if (
allocated(buf_i8_kind))
deallocate(buf_i8_kind)
300 if (
allocated(buf_r4_kind))
deallocate(buf_r4_kind)
301 if (
allocated(buf_r8_kind))
deallocate(buf_r8_kind)
311 corner, edge_lengths)
313 class(fmsnetcdffile_t),
intent(in) :: fileobj
314 character(len=*),
intent(in) :: variable_name
315 class(*),
contiguous,
intent(in),
target :: cdata(:,:,:)
319 integer,
intent(in),
optional :: unlim_dim_level
321 integer,
dimension(3),
intent(in),
optional :: corner
325 integer,
dimension(3),
intent(in),
optional :: edge_lengths
330 integer,
dimension(2) :: compressed_dim_index
333 integer,
dimension(3) :: c
334 integer,
dimension(3) :: e
336 integer(kind=i4_kind),
dimension(:,:,:),
allocatable :: buf_i4_kind
337 integer(kind=i8_kind),
dimension(:,:,:),
allocatable :: buf_i8_kind
338 real (kind=r4_kind),
dimension(:,:,:),
allocatable :: buf_r4_kind
339 real (kind=r8_kind),
dimension(:,:,:),
allocatable :: buf_r8_kind
341 character(len=200) :: append_error_msg
342 class(*),
pointer :: cdata_dummy(:,:,:,:)
351 append_error_msg =
"compressed_write_3d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
353 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
354 if (compressed_dim_index(1) .eq. dimension_not_found)
then
355 call netcdf_write_data(fileobj, variable_name, cdata, &
356 unlim_dim_level=unlim_dim_level, corner=corner, &
357 edge_lengths=edge_lengths)
359 else if (compressed_dim_index(1) .eq. 3)
then
360 cdata_dummy(1:
size(cdata,1), 1:
size(cdata,2), 1:
size(cdata,3), 1:1) => cdata(:,:,:)
366 if (fileobj%is_root)
then
367 e(compressed_dim_index(1)) = sum(fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems)
369 type is (
integer(kind=i4_kind))
370 call allocate_array(buf_i4_kind, e)
371 type is (
integer(kind=i8_kind))
372 call allocate_array(buf_i8_kind, e)
373 type is (real(kind=r4_kind))
374 call allocate_array(buf_r4_kind, e)
375 type is (real(kind=r8_kind))
376 call allocate_array(buf_r8_kind, e)
378 call error(
"unsupported variable type: "//trim(append_error_msg))
383 index = findloc(fileobj%pelist,
mpp_pe())
385 c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(index(1))
386 e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(index(1))
388 if (compressed_dim_index(1) .eq. 1)
then
389 is = c(compressed_dim_index(1))
390 ie = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
394 js = c(compressed_dim_index(1))
395 je = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
401 type is (
integer(kind=i4_kind))
402 call mpp_gather(is, ie, js, je,
size(cdata,3), &
403 fileobj%pelist, cdata, buf_i4_kind, fileobj%is_root)
404 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
405 unlim_dim_level=unlim_dim_level)
406 type is (
integer(kind=i8_kind))
407 call mpp_gather(is, ie, js, je,
size(cdata,3), &
408 fileobj%pelist, cdata, buf_i8_kind, fileobj%is_root)
409 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
410 unlim_dim_level=unlim_dim_level)
411 type is (real(kind=r4_kind))
412 call mpp_gather(is, ie, js, je,
size(cdata,3), &
413 fileobj%pelist, cdata, buf_r4_kind, fileobj%is_root)
414 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
415 unlim_dim_level=unlim_dim_level)
416 type is (real(kind=r8_kind))
417 call mpp_gather(is, ie, js, je,
size(cdata,3), &
418 fileobj%pelist, cdata, buf_r8_kind, fileobj%is_root)
419 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
420 unlim_dim_level=unlim_dim_level)
422 call error(
"unsupported variable type: "//trim(append_error_msg))
425 if (fileobj%is_root)
then
426 if (
allocated(buf_i4_kind))
deallocate(buf_i4_kind)
427 if (
allocated(buf_i8_kind))
deallocate(buf_i8_kind)
428 if (
allocated(buf_r4_kind))
deallocate(buf_r4_kind)
429 if (
allocated(buf_r8_kind))
deallocate(buf_r8_kind)
439 corner, edge_lengths)
441 class(fmsnetcdffile_t),
intent(in) :: fileobj
442 character(len=*),
intent(in) :: variable_name
443 class(*),
dimension(:,:,:,:),
intent(in) :: cdata
447 integer,
intent(in),
optional :: unlim_dim_level
449 integer,
dimension(4),
intent(in),
optional :: corner
453 integer,
dimension(4),
intent(in),
optional :: edge_lengths
458 integer,
dimension(2) :: compressed_dim_index
459 integer,
dimension(4) :: c
460 integer,
dimension(4) :: e
462 integer(kind=i4_kind),
dimension(:,:,:,:),
allocatable :: buf_i4_kind
463 integer(kind=i8_kind),
dimension(:,:,:,:),
allocatable :: buf_i8_kind
464 real(kind=r4_kind),
dimension(:,:,:,:),
allocatable :: buf_r4_kind
465 real(kind=r8_kind),
dimension(:,:,:,:),
allocatable :: buf_r8_kind
466 character(len=200) :: append_error_msg
468 append_error_msg =
"compressed_write_4d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
470 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
471 if (compressed_dim_index(1) .eq. dimension_not_found)
then
472 call netcdf_write_data(fileobj, variable_name, cdata, &
473 unlim_dim_level=unlim_dim_level, corner=corner, &
474 edge_lengths=edge_lengths)
479 if (fileobj%is_root)
then
482 do i = 1,
size(fileobj%pelist)
483 c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(i)
484 e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(i)
486 call netcdf_write_data(fileobj, variable_name, cdata, &
487 unlim_dim_level=unlim_dim_level, corner=c, &
491 type is (
integer(kind=i4_kind))
492 call allocate_array(buf_i4_kind, e)
493 call mpp_recv(buf_i4_kind,
size(buf_i4_kind), fileobj%pelist(i), block=.true.)
494 call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
495 unlim_dim_level=unlim_dim_level, corner=c, &
497 deallocate(buf_i4_kind)
498 type is (
integer(kind=i8_kind))
499 call allocate_array(buf_i8_kind, e)
500 call mpp_recv(buf_i8_kind,
size(buf_i8_kind), fileobj%pelist(i), block=.true.)
501 call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
502 unlim_dim_level=unlim_dim_level, corner=c, &
504 deallocate(buf_i8_kind)
505 type is (real(kind=r4_kind))
506 call allocate_array(buf_r4_kind, e)
507 call mpp_recv(buf_r4_kind,
size(buf_r4_kind), fileobj%pelist(i), block=.true.)
508 call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
509 unlim_dim_level=unlim_dim_level, corner=c, &
511 deallocate(buf_r4_kind)
512 type is (real(kind=r8_kind))
513 call allocate_array(buf_r8_kind, e)
514 call mpp_recv(buf_r8_kind,
size(buf_r8_kind), fileobj%pelist(i), block=.true.)
515 call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
516 unlim_dim_level=unlim_dim_level, corner=c, &
518 deallocate(buf_r8_kind)
520 call error(
"unsupported variable type: "//trim(append_error_msg))
526 type is (
integer(kind=i4_kind))
527 call mpp_send(cdata,
size(cdata), fileobj%io_root)
528 type is (
integer(kind=i8_kind))
529 call mpp_send(cdata,
size(cdata), fileobj%io_root)
530 type is (real(kind=r4_kind))
531 call mpp_send(cdata,
size(cdata), fileobj%io_root)
532 type is (real(kind=r8_kind))
533 call mpp_send(cdata,
size(cdata), fileobj%io_root)
535 call error(
"unsupported variable type: "//trim(append_error_msg))
547 corner, edge_lengths)
549 class(fmsnetcdffile_t),
intent(in) :: fileobj
550 character(len=*),
intent(in) :: variable_name
551 class(*),
dimension(:,:,:,:,:),
intent(in) :: cdata
555 integer,
intent(in),
optional :: unlim_dim_level
557 integer,
dimension(5),
intent(in),
optional :: corner
561 integer,
dimension(5),
intent(in),
optional :: edge_lengths
566 integer,
dimension(2) :: compressed_dim_index
567 integer,
dimension(5) :: c
568 integer,
dimension(5) :: e
570 integer(kind=i4_kind),
dimension(:,:,:,:,:),
allocatable :: buf_i4_kind
571 integer(kind=i8_kind),
dimension(:,:,:,:,:),
allocatable :: buf_i8_kind
572 real(kind=r4_kind),
dimension(:,:,:,:,:),
allocatable :: buf_r4_kind
573 real(kind=r8_kind),
dimension(:,:,:,:,:),
allocatable :: buf_r8_kind
574 character(len=200) :: append_error_msg
576 append_error_msg =
"compressed_write_5d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
578 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
579 if (compressed_dim_index(1) .eq. dimension_not_found)
then
580 call netcdf_write_data(fileobj, variable_name, cdata, &
581 unlim_dim_level=unlim_dim_level, corner=corner, &
582 edge_lengths=edge_lengths)
587 if (fileobj%is_root)
then
590 do i = 1,
size(fileobj%pelist)
591 c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(i)
592 e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(i)
594 call netcdf_write_data(fileobj, variable_name, cdata, &
595 unlim_dim_level=unlim_dim_level, corner=c, &
599 type is (
integer(kind=i4_kind))
600 call allocate_array(buf_i4_kind, e)
601 call mpp_recv(buf_i4_kind,
size(buf_i4_kind), fileobj%pelist(i), block=.true.)
602 call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
603 unlim_dim_level=unlim_dim_level, corner=c, &
605 deallocate(buf_i4_kind)
606 type is (
integer(kind=i8_kind))
607 call allocate_array(buf_i8_kind, e)
608 call mpp_recv(buf_i8_kind,
size(buf_i8_kind), fileobj%pelist(i), block=.true.)
609 call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
610 unlim_dim_level=unlim_dim_level, corner=c, &
612 deallocate(buf_i8_kind)
613 type is (real(kind=r4_kind))
614 call allocate_array(buf_r4_kind, e)
615 call mpp_recv(buf_r4_kind,
size(buf_r4_kind), fileobj%pelist(i), block=.true.)
616 call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
617 unlim_dim_level=unlim_dim_level, corner=c, &
619 deallocate(buf_r4_kind)
620 type is (real(kind=r8_kind))
621 call allocate_array(buf_r8_kind, e)
622 call mpp_recv(buf_r8_kind,
size(buf_r8_kind), fileobj%pelist(i), block=.true.)
623 call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
624 unlim_dim_level=unlim_dim_level, corner=c, &
626 deallocate(buf_r8_kind)
628 call error(
"unsupported variable type: "//trim(append_error_msg))
634 type is (
integer(kind=i4_kind))
635 call mpp_send(cdata,
size(cdata), fileobj%io_root)
636 type is (
integer(kind=i8_kind))
637 call mpp_send(cdata,
size(cdata), fileobj%io_root)
638 type is (real(kind=r4_kind))
639 call mpp_send(cdata,
size(cdata), fileobj%io_root)
640 type is (real(kind=r8_kind))
641 call mpp_send(cdata,
size(cdata), fileobj%io_root)
643 call error(
"unsupported variable type: "//trim(append_error_msg))
652 corner, edge_lengths)
654 type(fmsnetcdffile_t),
intent(in) :: fileobj
655 character(len=*),
intent(in) :: variable_name
656 class(*),
dimension(:),
intent(in) :: cdata
660 integer,
intent(in),
optional :: unlim_dim_level
662 integer,
dimension(1),
intent(in),
optional :: corner
666 integer,
dimension(1),
intent(in),
optional :: edge_lengths
671 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
672 edge_lengths=edge_lengths)
679 corner, edge_lengths)
681 type(fmsnetcdffile_t),
intent(in) :: fileobj
682 character(len=*),
intent(in) :: variable_name
683 class(*),
dimension(:,:),
intent(in) :: cdata
687 integer,
intent(in),
optional :: unlim_dim_level
689 integer,
dimension(2),
intent(in),
optional :: corner
693 integer,
dimension(2),
intent(in),
optional :: edge_lengths
698 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
699 edge_lengths=edge_lengths)
705 corner, edge_lengths)
707 type(fmsnetcdffile_t),
intent(in) :: fileobj
708 character(len=*),
intent(in) :: variable_name
709 class(*),
dimension(:,:,:),
intent(in) :: cdata
713 integer,
intent(in),
optional :: unlim_dim_level
715 integer,
dimension(3),
intent(in),
optional :: corner
719 integer,
dimension(3),
intent(in),
optional :: edge_lengths
724 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
725 edge_lengths=edge_lengths)
731 corner, edge_lengths)
733 type(fmsnetcdffile_t),
intent(in) :: fileobj
734 character(len=*),
intent(in) :: variable_name
735 class(*),
dimension(:,:,:,:),
intent(in) :: cdata
739 integer,
intent(in),
optional :: unlim_dim_level
741 integer,
dimension(4),
intent(in),
optional :: corner
745 integer,
dimension(4),
intent(in),
optional :: edge_lengths
750 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
751 edge_lengths=edge_lengths)
757 corner, edge_lengths)
759 type(fmsnetcdffile_t),
intent(in) :: fileobj
760 character(len=*),
intent(in) :: variable_name
761 class(*),
dimension(:,:,:,:,:),
intent(in) :: cdata
765 integer,
intent(in),
optional :: unlim_dim_level
767 integer,
dimension(5),
intent(in),
optional :: corner
771 integer,
dimension(5),
intent(in),
optional :: edge_lengths
776 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
777 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.