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))
144 type is (
integer(kind=i4_kind))
145 call mpp_gather(cdata,
size(cdata), buf_i4_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
147 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
148 unlim_dim_level=unlim_dim_level)
149 type is (
integer(kind=i8_kind))
150 call mpp_gather(cdata,
size(cdata), buf_i8_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
152 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
153 unlim_dim_level=unlim_dim_level)
154 type is (real(kind=r4_kind))
155 call mpp_gather(cdata,
size(cdata), buf_r4_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
157 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
158 unlim_dim_level=unlim_dim_level)
159 type is (real(kind=r8_kind))
160 call mpp_gather(cdata,
size(cdata), buf_r8_kind, fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems, &
162 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
163 unlim_dim_level=unlim_dim_level)
165 call error(
"unsupported variable type: "//trim(append_error_msg))
167 if (fileobj%is_root)
then
168 if (
allocated(buf_i4_kind))
deallocate(buf_i4_kind)
169 if (
allocated(buf_i8_kind))
deallocate(buf_i8_kind)
170 if (
allocated(buf_r4_kind))
deallocate(buf_r4_kind)
171 if (
allocated(buf_r8_kind))
deallocate(buf_r8_kind)
181 corner, edge_lengths)
183 class(fmsnetcdffile_t),
intent(in) :: fileobj
184 character(len=*),
intent(in) :: variable_name
185 class(*),
dimension(:,:),
intent(in) :: cdata
189 integer,
intent(in),
optional :: unlim_dim_level
191 integer,
dimension(2),
intent(in),
optional :: corner
195 integer,
dimension(2),
intent(in),
optional :: edge_lengths
200 integer,
dimension(2) :: compressed_dim_index
203 integer,
dimension(2) :: c
204 integer,
dimension(2) :: e
206 integer(kind=i4_kind),
dimension(:,:),
allocatable :: buf_i4_kind
207 integer(kind=i8_kind),
dimension(:,:),
allocatable :: buf_i8_kind
208 real (kind=r4_kind),
dimension(:,:),
allocatable :: buf_r4_kind
209 real (kind=r8_kind),
dimension(:,:),
allocatable :: buf_r8_kind
211 character(len=200) :: append_error_msg
220 append_error_msg =
"compressed_write_2d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
222 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
223 if (compressed_dim_index(1) .eq. dimension_not_found)
then
224 call netcdf_write_data(fileobj, variable_name, cdata, &
225 unlim_dim_level=unlim_dim_level, corner=corner, &
226 edge_lengths=edge_lengths)
232 if (fileobj%is_root)
then
233 e(compressed_dim_index(1)) = sum(fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems)
235 type is (
integer(kind=i4_kind))
236 call allocate_array(buf_i4_kind, e)
237 type is (
integer(kind=i8_kind))
238 call allocate_array(buf_i8_kind, e)
239 type is (real(kind=r4_kind))
240 call allocate_array(buf_r4_kind, e)
241 type is (real(kind=r8_kind))
242 call allocate_array(buf_r8_kind, e)
244 call error(
"unsupported variable type: "//trim(append_error_msg))
249 index = findloc(fileobj%pelist,
mpp_pe())
251 c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(index(1))
252 e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(index(1))
254 if (compressed_dim_index(1) .eq. 1)
then
255 is = c(compressed_dim_index(1))
256 ie = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
260 js = c(compressed_dim_index(1))
261 je = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
267 type is (
integer(kind=i4_kind))
268 call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_i4_kind, fileobj%is_root)
269 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
270 unlim_dim_level=unlim_dim_level)
271 type is (
integer(kind=i8_kind))
272 call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_i8_kind, fileobj%is_root)
273 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
274 unlim_dim_level=unlim_dim_level)
275 type is (real(kind=r4_kind))
276 call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_r4_kind, fileobj%is_root)
277 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
278 unlim_dim_level=unlim_dim_level)
279 type is (real(kind=r8_kind))
280 call mpp_gather(is, ie, js, je, fileobj%pelist, cdata, buf_r8_kind, fileobj%is_root)
281 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
282 unlim_dim_level=unlim_dim_level)
284 call error(
"unsupported variable type: "//trim(append_error_msg))
287 if (fileobj%is_root)
then
288 if (
allocated(buf_i4_kind))
deallocate(buf_i4_kind)
289 if (
allocated(buf_i8_kind))
deallocate(buf_i8_kind)
290 if (
allocated(buf_r4_kind))
deallocate(buf_r4_kind)
291 if (
allocated(buf_r8_kind))
deallocate(buf_r8_kind)
301 corner, edge_lengths)
303 class(fmsnetcdffile_t),
intent(in) :: fileobj
304 character(len=*),
intent(in) :: variable_name
305 class(*),
contiguous,
intent(in),
target :: cdata(:,:,:)
309 integer,
intent(in),
optional :: unlim_dim_level
311 integer,
dimension(3),
intent(in),
optional :: corner
315 integer,
dimension(3),
intent(in),
optional :: edge_lengths
320 integer,
dimension(2) :: compressed_dim_index
323 integer,
dimension(3) :: c
324 integer,
dimension(3) :: e
326 integer(kind=i4_kind),
dimension(:,:,:),
allocatable :: buf_i4_kind
327 integer(kind=i8_kind),
dimension(:,:,:),
allocatable :: buf_i8_kind
328 real (kind=r4_kind),
dimension(:,:,:),
allocatable :: buf_r4_kind
329 real (kind=r8_kind),
dimension(:,:,:),
allocatable :: buf_r8_kind
331 character(len=200) :: append_error_msg
332 class(*),
pointer :: cdata_dummy(:,:,:,:)
341 append_error_msg =
"compressed_write_3d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
343 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
344 if (compressed_dim_index(1) .eq. dimension_not_found)
then
345 call netcdf_write_data(fileobj, variable_name, cdata, &
346 unlim_dim_level=unlim_dim_level, corner=corner, &
347 edge_lengths=edge_lengths)
349 else if (compressed_dim_index(1) .eq. 3)
then
350 cdata_dummy(1:
size(cdata,1), 1:
size(cdata,2), 1:
size(cdata,3), 1:1) => cdata(:,:,:)
356 if (fileobj%is_root)
then
357 e(compressed_dim_index(1)) = sum(fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems)
359 type is (
integer(kind=i4_kind))
360 call allocate_array(buf_i4_kind, e)
361 type is (
integer(kind=i8_kind))
362 call allocate_array(buf_i8_kind, e)
363 type is (real(kind=r4_kind))
364 call allocate_array(buf_r4_kind, e)
365 type is (real(kind=r8_kind))
366 call allocate_array(buf_r8_kind, e)
368 call error(
"unsupported variable type: "//trim(append_error_msg))
373 index = findloc(fileobj%pelist,
mpp_pe())
375 c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(index(1))
376 e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(index(1))
378 if (compressed_dim_index(1) .eq. 1)
then
379 is = c(compressed_dim_index(1))
380 ie = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
384 js = c(compressed_dim_index(1))
385 je = c(compressed_dim_index(1)) + e(compressed_dim_index(1)) - 1
391 type is (
integer(kind=i4_kind))
392 call mpp_gather(is, ie, js, je,
size(cdata,3), &
393 fileobj%pelist, cdata, buf_i4_kind, fileobj%is_root)
394 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
395 unlim_dim_level=unlim_dim_level)
396 type is (
integer(kind=i8_kind))
397 call mpp_gather(is, ie, js, je,
size(cdata,3), &
398 fileobj%pelist, cdata, buf_i8_kind, fileobj%is_root)
399 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
400 unlim_dim_level=unlim_dim_level)
401 type is (real(kind=r4_kind))
402 call mpp_gather(is, ie, js, je,
size(cdata,3), &
403 fileobj%pelist, cdata, buf_r4_kind, fileobj%is_root)
404 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
405 unlim_dim_level=unlim_dim_level)
406 type is (real(kind=r8_kind))
407 call mpp_gather(is, ie, js, je,
size(cdata,3), &
408 fileobj%pelist, cdata, buf_r8_kind, fileobj%is_root)
409 if (fileobj%is_root)
call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
410 unlim_dim_level=unlim_dim_level)
412 call error(
"unsupported variable type: "//trim(append_error_msg))
415 if (fileobj%is_root)
then
416 if (
allocated(buf_i4_kind))
deallocate(buf_i4_kind)
417 if (
allocated(buf_i8_kind))
deallocate(buf_i8_kind)
418 if (
allocated(buf_r4_kind))
deallocate(buf_r4_kind)
419 if (
allocated(buf_r8_kind))
deallocate(buf_r8_kind)
429 corner, edge_lengths)
431 class(fmsnetcdffile_t),
intent(in) :: fileobj
432 character(len=*),
intent(in) :: variable_name
433 class(*),
dimension(:,:,:,:),
intent(in) :: cdata
437 integer,
intent(in),
optional :: unlim_dim_level
439 integer,
dimension(4),
intent(in),
optional :: corner
443 integer,
dimension(4),
intent(in),
optional :: edge_lengths
448 integer,
dimension(2) :: compressed_dim_index
449 integer,
dimension(4) :: c
450 integer,
dimension(4) :: e
452 integer(kind=i4_kind),
dimension(:,:,:,:),
allocatable :: buf_i4_kind
453 integer(kind=i8_kind),
dimension(:,:,:,:),
allocatable :: buf_i8_kind
454 real(kind=r4_kind),
dimension(:,:,:,:),
allocatable :: buf_r4_kind
455 real(kind=r8_kind),
dimension(:,:,:,:),
allocatable :: buf_r8_kind
456 character(len=200) :: append_error_msg
458 append_error_msg =
"compressed_write_4d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
460 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
461 if (compressed_dim_index(1) .eq. dimension_not_found)
then
462 call netcdf_write_data(fileobj, variable_name, cdata, &
463 unlim_dim_level=unlim_dim_level, corner=corner, &
464 edge_lengths=edge_lengths)
469 if (fileobj%is_root)
then
472 do i = 1,
size(fileobj%pelist)
473 c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(i)
474 e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(i)
476 call netcdf_write_data(fileobj, variable_name, cdata, &
477 unlim_dim_level=unlim_dim_level, corner=c, &
481 type is (
integer(kind=i4_kind))
482 call allocate_array(buf_i4_kind, e)
483 call mpp_recv(buf_i4_kind,
size(buf_i4_kind), fileobj%pelist(i), block=.true.)
484 call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
485 unlim_dim_level=unlim_dim_level, corner=c, &
487 deallocate(buf_i4_kind)
488 type is (
integer(kind=i8_kind))
489 call allocate_array(buf_i8_kind, e)
490 call mpp_recv(buf_i8_kind,
size(buf_i8_kind), fileobj%pelist(i), block=.true.)
491 call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
492 unlim_dim_level=unlim_dim_level, corner=c, &
494 deallocate(buf_i8_kind)
495 type is (real(kind=r4_kind))
496 call allocate_array(buf_r4_kind, e)
497 call mpp_recv(buf_r4_kind,
size(buf_r4_kind), fileobj%pelist(i), block=.true.)
498 call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
499 unlim_dim_level=unlim_dim_level, corner=c, &
501 deallocate(buf_r4_kind)
502 type is (real(kind=r8_kind))
503 call allocate_array(buf_r8_kind, e)
504 call mpp_recv(buf_r8_kind,
size(buf_r8_kind), fileobj%pelist(i), block=.true.)
505 call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
506 unlim_dim_level=unlim_dim_level, corner=c, &
508 deallocate(buf_r8_kind)
510 call error(
"unsupported variable type: "//trim(append_error_msg))
516 type is (
integer(kind=i4_kind))
517 call mpp_send(cdata,
size(cdata), fileobj%io_root)
518 type is (
integer(kind=i8_kind))
519 call mpp_send(cdata,
size(cdata), fileobj%io_root)
520 type is (real(kind=r4_kind))
521 call mpp_send(cdata,
size(cdata), fileobj%io_root)
522 type is (real(kind=r8_kind))
523 call mpp_send(cdata,
size(cdata), fileobj%io_root)
525 call error(
"unsupported variable type: "//trim(append_error_msg))
537 corner, edge_lengths)
539 class(fmsnetcdffile_t),
intent(in) :: fileobj
540 character(len=*),
intent(in) :: variable_name
541 class(*),
dimension(:,:,:,:,:),
intent(in) :: cdata
545 integer,
intent(in),
optional :: unlim_dim_level
547 integer,
dimension(5),
intent(in),
optional :: corner
551 integer,
dimension(5),
intent(in),
optional :: edge_lengths
556 integer,
dimension(2) :: compressed_dim_index
557 integer,
dimension(5) :: c
558 integer,
dimension(5) :: e
560 integer(kind=i4_kind),
dimension(:,:,:,:,:),
allocatable :: buf_i4_kind
561 integer(kind=i8_kind),
dimension(:,:,:,:,:),
allocatable :: buf_i8_kind
562 real(kind=r4_kind),
dimension(:,:,:,:,:),
allocatable :: buf_r4_kind
563 real(kind=r8_kind),
dimension(:,:,:,:,:),
allocatable :: buf_r8_kind
564 character(len=200) :: append_error_msg
566 append_error_msg =
"compressed_write_5d: file:"//trim(fileobj%path)//
" variable:"//trim(variable_name)
568 compressed_dim_index = get_variable_compressed_dimension_index(fileobj, variable_name)
569 if (compressed_dim_index(1) .eq. dimension_not_found)
then
570 call netcdf_write_data(fileobj, variable_name, cdata, &
571 unlim_dim_level=unlim_dim_level, corner=corner, &
572 edge_lengths=edge_lengths)
577 if (fileobj%is_root)
then
580 do i = 1,
size(fileobj%pelist)
581 c(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_corner(i)
582 e(compressed_dim_index(1)) = fileobj%compressed_dims(compressed_dim_index(2))%npes_nelems(i)
584 call netcdf_write_data(fileobj, variable_name, cdata, &
585 unlim_dim_level=unlim_dim_level, corner=c, &
589 type is (
integer(kind=i4_kind))
590 call allocate_array(buf_i4_kind, e)
591 call mpp_recv(buf_i4_kind,
size(buf_i4_kind), fileobj%pelist(i), block=.true.)
592 call netcdf_write_data(fileobj, variable_name, buf_i4_kind, &
593 unlim_dim_level=unlim_dim_level, corner=c, &
595 deallocate(buf_i4_kind)
596 type is (
integer(kind=i8_kind))
597 call allocate_array(buf_i8_kind, e)
598 call mpp_recv(buf_i8_kind,
size(buf_i8_kind), fileobj%pelist(i), block=.true.)
599 call netcdf_write_data(fileobj, variable_name, buf_i8_kind, &
600 unlim_dim_level=unlim_dim_level, corner=c, &
602 deallocate(buf_i8_kind)
603 type is (real(kind=r4_kind))
604 call allocate_array(buf_r4_kind, e)
605 call mpp_recv(buf_r4_kind,
size(buf_r4_kind), fileobj%pelist(i), block=.true.)
606 call netcdf_write_data(fileobj, variable_name, buf_r4_kind, &
607 unlim_dim_level=unlim_dim_level, corner=c, &
609 deallocate(buf_r4_kind)
610 type is (real(kind=r8_kind))
611 call allocate_array(buf_r8_kind, e)
612 call mpp_recv(buf_r8_kind,
size(buf_r8_kind), fileobj%pelist(i), block=.true.)
613 call netcdf_write_data(fileobj, variable_name, buf_r8_kind, &
614 unlim_dim_level=unlim_dim_level, corner=c, &
616 deallocate(buf_r8_kind)
618 call error(
"unsupported variable type: "//trim(append_error_msg))
624 type is (
integer(kind=i4_kind))
625 call mpp_send(cdata,
size(cdata), fileobj%io_root)
626 type is (
integer(kind=i8_kind))
627 call mpp_send(cdata,
size(cdata), fileobj%io_root)
628 type is (real(kind=r4_kind))
629 call mpp_send(cdata,
size(cdata), fileobj%io_root)
630 type is (real(kind=r8_kind))
631 call mpp_send(cdata,
size(cdata), fileobj%io_root)
633 call error(
"unsupported variable type: "//trim(append_error_msg))
642 corner, edge_lengths)
644 type(fmsnetcdffile_t),
intent(in) :: fileobj
645 character(len=*),
intent(in) :: variable_name
646 class(*),
dimension(:),
intent(in) :: cdata
650 integer,
intent(in),
optional :: unlim_dim_level
652 integer,
dimension(1),
intent(in),
optional :: corner
656 integer,
dimension(1),
intent(in),
optional :: edge_lengths
661 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
662 edge_lengths=edge_lengths)
669 corner, edge_lengths)
671 type(fmsnetcdffile_t),
intent(in) :: fileobj
672 character(len=*),
intent(in) :: variable_name
673 class(*),
dimension(:,:),
intent(in) :: cdata
677 integer,
intent(in),
optional :: unlim_dim_level
679 integer,
dimension(2),
intent(in),
optional :: corner
683 integer,
dimension(2),
intent(in),
optional :: edge_lengths
688 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
689 edge_lengths=edge_lengths)
695 corner, edge_lengths)
697 type(fmsnetcdffile_t),
intent(in) :: fileobj
698 character(len=*),
intent(in) :: variable_name
699 class(*),
dimension(:,:,:),
intent(in) :: cdata
703 integer,
intent(in),
optional :: unlim_dim_level
705 integer,
dimension(3),
intent(in),
optional :: corner
709 integer,
dimension(3),
intent(in),
optional :: edge_lengths
714 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
715 edge_lengths=edge_lengths)
721 corner, edge_lengths)
723 type(fmsnetcdffile_t),
intent(in) :: fileobj
724 character(len=*),
intent(in) :: variable_name
725 class(*),
dimension(:,:,:,:),
intent(in) :: cdata
729 integer,
intent(in),
optional :: unlim_dim_level
731 integer,
dimension(4),
intent(in),
optional :: corner
735 integer,
dimension(4),
intent(in),
optional :: edge_lengths
740 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
741 edge_lengths=edge_lengths)
747 corner, edge_lengths)
749 type(fmsnetcdffile_t),
intent(in) :: fileobj
750 character(len=*),
intent(in) :: variable_name
751 class(*),
dimension(:,:,:,:,:),
intent(in) :: cdata
755 integer,
intent(in),
optional :: unlim_dim_level
757 integer,
dimension(5),
intent(in),
optional :: corner
761 integer,
dimension(5),
intent(in),
optional :: edge_lengths
766 call compressed_write(fileobj, variable_name, cdata, unlim_dim_level, corner, &
767 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.