29 subroutine domain_read_0d(fileobj, variable_name, vdata, unlim_dim_level, corner)
31 type(fmsnetcdfdomainfile_t),
intent(in) :: fileobj
32 character(len=*),
intent(in) :: variable_name
33 class(*),
intent(inout) :: vdata
36 integer,
intent(in),
optional :: unlim_dim_level
38 integer,
intent(in),
optional :: corner
43 call netcdf_read_data(fileobj, variable_name, vdata, &
44 unlim_dim_level=unlim_dim_level, corner=corner, &
58 type(fmsnetcdfdomainfile_t),
intent(in) :: fileobj
59 character(len=*),
intent(in) :: variable_name
60 class(*),
dimension(:),
intent(inout) :: vdata
63 integer,
intent(in),
optional :: unlim_dim_level
65 integer,
dimension(1),
intent(in),
optional :: corner
69 integer,
dimension(1),
intent(in),
optional :: edge_lengths
74 call netcdf_read_data(fileobj, variable_name, vdata, &
75 unlim_dim_level=unlim_dim_level, corner=corner, &
76 edge_lengths=edge_lengths, broadcast=.true.)
88 type(fmsnetcdfdomainfile_t),
intent(in) :: fileobj
89 character(len=*),
intent(in) :: variable_name
90 class(*),
contiguous,
target,
intent(inout) :: vdata(:,:)
93 integer,
intent(in),
optional :: unlim_dim_level
95 integer,
dimension(2),
intent(in),
optional :: corner
99 integer,
dimension(2),
intent(in),
optional :: edge_lengths
104 integer :: xdim_index
105 integer :: ydim_index
117 logical :: buffer_includes_halos
122 type(domain2d),
pointer :: io_domain
125 integer(kind=i4_kind),
dimension(:,:),
allocatable :: buf_i4_kind_pe
126 integer(kind=i8_kind),
dimension(:,:),
allocatable :: buf_i8_kind_pe
127 real(kind=r4_kind),
dimension(:,:),
allocatable :: buf_r4_kind_pe
128 real(kind=r8_kind),
dimension(:,:),
allocatable :: buf_r8_kind_pe
129 integer(kind=i4_kind),
dimension(:,:),
allocatable :: buf_i4_kind
130 integer(kind=i8_kind),
dimension(:,:),
allocatable :: buf_i8_kind
131 real(kind=r4_kind),
dimension(:,:),
allocatable :: buf_r4_kind
132 real(kind=r8_kind),
dimension(:,:),
allocatable :: buf_r8_kind
133 class(*),
dimension(:,:,:,:),
pointer :: vdata_dummy
135 if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
136 xdim_index, ydim_index, xpos, ypos))
then
137 call netcdf_read_data(fileobj, variable_name, vdata, &
138 unlim_dim_level=unlim_dim_level, corner=corner, &
139 edge_lengths=edge_lengths, broadcast=.true.)
143 if (xdim_index .ne. 1 .or. ydim_index .ne. 2)
then
147 vdata_dummy(1:
size(vdata,1),1:
size(vdata,2), 1:1, 1:1) => vdata(:,:)
148 call domain_read_4d(fileobj, variable_name, vdata_dummy, unlim_dim_level)
152 call domain_offsets(
size(vdata, xdim_index),
size(vdata, ydim_index), fileobj%domain, &
153 xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
154 msg=
"file:"//trim(fileobj%path)//
" and variable:"//trim(variable_name))
158 call mpp_get_global_domain(io_domain, xbegin=xgbegin, xsize=xgsize, position=xpos)
159 call mpp_get_global_domain(io_domain, ybegin=ygbegin, ysize=ygsize, position=ypos)
162 if (fileobj%is_root)
then
164 if (fileobj%adjust_indices)
then
171 c(xdim_index) = xgbegin
172 c(ydim_index) = ygbegin
175 e(xdim_index) = xgsize
176 e(ydim_index) = ygsize
180 type is (
integer(kind=i4_kind))
181 call allocate_array(buf_i4_kind, e)
182 call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
183 unlim_dim_level=unlim_dim_level, &
184 corner=c, edge_lengths=e, broadcast=.false.)
185 type is (
integer(kind=i8_kind))
186 call allocate_array(buf_i8_kind, e)
187 call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
188 unlim_dim_level=unlim_dim_level, &
189 corner=c, edge_lengths=e, broadcast=.false.)
190 type is (real(kind=r4_kind))
191 call allocate_array(buf_r4_kind, e)
192 call netcdf_read_data(fileobj, variable_name, buf_r4_kind, &
193 unlim_dim_level=unlim_dim_level, &
194 corner=c, edge_lengths=e, broadcast=.false.)
195 type is (real(kind=r8_kind))
196 call allocate_array(buf_r8_kind, e)
197 call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
198 unlim_dim_level=unlim_dim_level, &
199 corner=c, edge_lengths=e, broadcast=.false.)
201 call error(
"unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//
" variable:"// &
202 & trim(variable_name))
210 if (buffer_includes_halos)
then
212 c(xdim_index) = isc - isd + 1
213 c(ydim_index) = jsc - jsd + 1
219 e(xdim_index) = xc_size
220 e(ydim_index) = yc_size
223 type is (
integer(kind=i4_kind))
224 call allocate_array(buf_i4_kind_pe, e)
225 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, fileobj%pelist, &
226 buf_i4_kind_pe, buf_i4_kind, fileobj%is_root)
227 call put_array_section(buf_i4_kind_pe, vdata, c, e)
228 deallocate(buf_i4_kind_pe)
229 type is (
integer(kind=i8_kind))
230 call allocate_array(buf_i8_kind_pe, e)
231 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, fileobj%pelist, &
232 buf_i8_kind_pe, buf_i8_kind, fileobj%is_root)
233 call put_array_section(buf_i8_kind_pe, vdata, c, e)
234 deallocate(buf_i8_kind_pe)
235 type is (real(kind=r4_kind))
236 call allocate_array(buf_r4_kind_pe, e)
237 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, fileobj%pelist, &
238 buf_r4_kind_pe, buf_r4_kind, fileobj%is_root)
239 call put_array_section(buf_r4_kind_pe, vdata, c, e)
240 deallocate(buf_r4_kind_pe)
241 type is (real(kind=r8_kind))
242 call allocate_array(buf_r8_kind_pe, e)
243 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, fileobj%pelist, &
244 buf_r8_kind_pe, buf_r8_kind, fileobj%is_root)
245 call put_array_section(buf_r8_kind_pe, vdata, c, e)
246 deallocate(buf_r8_kind_pe)
248 call error(
"unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//
" variable:"// &
249 & trim(variable_name))
252 if (fileobj%is_root)
then
253 if (
allocated(buf_i4_kind))
deallocate(buf_i4_kind)
254 if (
allocated(buf_i8_kind))
deallocate(buf_i8_kind)
255 if (
allocated(buf_r4_kind))
deallocate(buf_r4_kind)
256 if (
allocated(buf_r8_kind))
deallocate(buf_r8_kind)
267 corner, edge_lengths)
268 type(fmsnetcdfdomainfile_t),
intent(in) :: fileobj
269 character(len=*),
intent(in) :: variable_name
270 class(*),
contiguous,
target,
intent(inout) :: vdata(:,:,:)
273 integer,
intent(in),
optional :: unlim_dim_level
275 integer,
dimension(3),
intent(in),
optional :: corner
279 integer,
dimension(3),
intent(in),
optional :: edge_lengths
284 integer :: xdim_index
285 integer :: ydim_index
297 logical :: buffer_includes_halos
302 type(domain2d),
pointer :: io_domain
305 integer(kind=i4_kind),
dimension(:,:,:),
allocatable :: buf_i4_kind_pe
306 integer(kind=i8_kind),
dimension(:,:,:),
allocatable :: buf_i8_kind_pe
307 real(kind=r4_kind),
dimension(:,:,:),
allocatable :: buf_r4_kind_pe
308 real(kind=r8_kind),
dimension(:,:,:),
allocatable :: buf_r8_kind_pe
309 integer(kind=i4_kind),
dimension(:,:,:),
allocatable :: buf_i4_kind
310 integer(kind=i8_kind),
dimension(:,:,:),
allocatable :: buf_i8_kind
311 real(kind=r4_kind),
dimension(:,:,:),
allocatable :: buf_r4_kind
312 real(kind=r8_kind),
dimension(:,:,:),
allocatable :: buf_r8_kind
313 class(*),
dimension(:,:,:,:),
pointer :: vdata_dummy
315 if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
316 xdim_index, ydim_index, xpos, ypos))
then
317 call netcdf_read_data(fileobj, variable_name, vdata, &
318 unlim_dim_level=unlim_dim_level, corner=corner, &
319 edge_lengths=edge_lengths, broadcast=.true.)
323 if (xdim_index .ne. 1 .or. ydim_index .ne. 2)
then
327 vdata_dummy(1:
size(vdata,1),1:
size(vdata,2), 1:
size(vdata,3), 1:1) => vdata(:,:,:)
328 call domain_read_4d(fileobj, variable_name, vdata_dummy, unlim_dim_level)
332 call domain_offsets(
size(vdata, xdim_index),
size(vdata, ydim_index), fileobj%domain, &
333 xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
334 msg=
"file:"//trim(fileobj%path)//
" and variable:"//trim(variable_name))
336 if (
present(corner)) c = corner
339 if (
present(edge_lengths)) e = edge_lengths
341 call mpp_get_global_domain(io_domain, xbegin=xgbegin, xsize=xgsize, position=xpos)
342 call mpp_get_global_domain(io_domain, ybegin=ygbegin, ysize=ygsize, position=ypos)
345 if (fileobj%is_root)
then
347 if (fileobj%adjust_indices)
then
354 c(xdim_index) = xgbegin
355 c(ydim_index) = ygbegin
358 e(xdim_index) = xgsize
359 e(ydim_index) = ygsize
363 type is (
integer(kind=i4_kind))
364 call allocate_array(buf_i4_kind, e)
365 call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
366 unlim_dim_level=unlim_dim_level, &
367 corner=c, edge_lengths=e, broadcast=.false.)
368 type is (
integer(kind=i8_kind))
369 call allocate_array(buf_i8_kind, e)
370 call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
371 unlim_dim_level=unlim_dim_level, &
372 corner=c, edge_lengths=e, broadcast=.false.)
373 type is (real(kind=r4_kind))
374 call allocate_array(buf_r4_kind, e)
375 call netcdf_read_data(fileobj, variable_name, buf_r4_kind, &
376 unlim_dim_level=unlim_dim_level, &
377 corner=c, edge_lengths=e, broadcast=.false.)
378 type is (real(kind=r8_kind))
379 call allocate_array(buf_r8_kind, e)
380 call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
381 unlim_dim_level=unlim_dim_level, &
382 corner=c, edge_lengths=e, broadcast=.false.)
384 call error(
"unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//
" variable:"// &
385 & trim(variable_name))
393 if (buffer_includes_halos)
then
395 c(xdim_index) = isc - isd + 1
396 c(ydim_index) = jsc - jsd + 1
402 e(xdim_index) = xc_size
403 e(ydim_index) = yc_size
406 type is (
integer(kind=i4_kind))
407 call allocate_array(buf_i4_kind_pe, e)
408 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
409 buf_i4_kind_pe, buf_i4_kind, fileobj%is_root)
410 call put_array_section(buf_i4_kind_pe, vdata, c, e)
411 deallocate(buf_i4_kind_pe)
412 type is (
integer(kind=i8_kind))
413 call allocate_array(buf_i8_kind_pe, e)
414 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
415 buf_i8_kind_pe, buf_i8_kind, fileobj%is_root)
416 call put_array_section(buf_i8_kind_pe, vdata, c, e)
417 deallocate(buf_i8_kind_pe)
418 type is (real(kind=r4_kind))
419 call allocate_array(buf_r4_kind_pe, e)
420 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
421 buf_r4_kind_pe, buf_r4_kind, fileobj%is_root)
422 call put_array_section(buf_r4_kind_pe, vdata, c, e)
423 deallocate(buf_r4_kind_pe)
424 type is (real(kind=r8_kind))
425 call allocate_array(buf_r8_kind_pe, e)
426 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
427 buf_r8_kind_pe, buf_r8_kind, fileobj%is_root)
428 call put_array_section(buf_r8_kind_pe, vdata, c, e)
429 deallocate(buf_r8_kind_pe)
431 call error(
"unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//
" variable:"// &
432 & trim(variable_name))
435 if (fileobj%is_root)
then
436 if (
allocated(buf_i4_kind))
deallocate(buf_i4_kind)
437 if (
allocated(buf_i8_kind))
deallocate(buf_i8_kind)
438 if (
allocated(buf_r4_kind))
deallocate(buf_r4_kind)
439 if (
allocated(buf_r8_kind))
deallocate(buf_r8_kind)
451 corner, edge_lengths)
453 type(fmsnetcdfdomainfile_t),
intent(in) :: fileobj
454 character(len=*),
intent(in) :: variable_name
455 class(*),
dimension(:,:,:,:),
intent(inout) :: vdata
458 integer,
intent(in),
optional :: unlim_dim_level
460 integer,
dimension(4),
intent(in),
optional :: corner
464 integer,
dimension(4),
intent(in),
optional :: edge_lengths
469 integer :: xdim_index
470 integer :: ydim_index
471 type(domain2d),
pointer :: io_domain
481 integer,
dimension(:),
allocatable :: pe_isc
482 integer,
dimension(:),
allocatable :: pe_icsize
483 integer,
dimension(:),
allocatable :: pe_jsc
484 integer,
dimension(:),
allocatable :: pe_jcsize
485 integer,
dimension(4) :: c
486 integer,
dimension(4) :: e
487 integer(kind=i4_kind),
dimension(:,:,:,:),
allocatable :: buf_i4_kind
488 integer(kind=i8_kind),
dimension(:,:,:,:),
allocatable :: buf_i8_kind
489 real(kind=r4_kind),
dimension(:,:,:,:),
allocatable :: buf_r4_kind
490 real(kind=r8_kind),
dimension(:,:,:,:),
allocatable :: buf_r8_kind
491 logical :: buffer_includes_halos
495 if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
496 xdim_index, ydim_index, xpos, ypos))
then
497 call netcdf_read_data(fileobj, variable_name, vdata, &
498 unlim_dim_level=unlim_dim_level, corner=corner, &
499 edge_lengths=edge_lengths, broadcast=.true.)
503 call domain_offsets(
size(vdata, xdim_index),
size(vdata, ydim_index), fileobj%domain, &
504 xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
505 msg=
"file:"//trim(fileobj%path)//
" and variable:"//trim(variable_name))
508 if (
present(edge_lengths)) e = edge_lengths
511 if (fileobj%is_root)
then
512 allocate(pe_isc(
size(fileobj%pelist)))
513 allocate(pe_icsize(
size(fileobj%pelist)))
514 allocate(pe_jsc(
size(fileobj%pelist)))
515 allocate(pe_jcsize(
size(fileobj%pelist)))
516 call mpp_get_compute_domains(io_domain, xbegin=pe_isc, xsize=pe_icsize, position=xpos)
517 call mpp_get_compute_domains(io_domain, ybegin=pe_jsc, ysize=pe_jcsize, position=ypos)
518 call mpp_get_global_domain(io_domain, xbegin=xgmin, position=xpos)
519 call mpp_get_global_domain(io_domain, ybegin=ygmin, position=ypos)
520 do i = 1,
size(fileobj%pelist)
521 if (
present(corner)) c = corner
522 c(xdim_index) = pe_isc(i)
523 c(ydim_index) = pe_jsc(i)
524 if (fileobj%adjust_indices)
then
525 c(xdim_index) = c(xdim_index) - xgmin + 1
526 c(ydim_index) = c(ydim_index) - ygmin + 1
528 e(xdim_index) = pe_icsize(i)
529 e(ydim_index) = pe_jcsize(i)
531 type is (
integer(kind=i4_kind))
533 call allocate_array(buf_i4_kind, e)
534 call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
535 unlim_dim_level=unlim_dim_level, &
536 corner=c, edge_lengths=e, broadcast=.false.)
540 if (buffer_includes_halos)
then
542 c(xdim_index) = isc - isd + 1
543 c(ydim_index) = jsc - jsd + 1
545 call put_array_section(buf_i4_kind, vdata, c, e)
548 call mpp_send(buf_i4_kind,
size(buf_i4_kind), fileobj%pelist(i))
551 deallocate(buf_i4_kind)
552 type is (
integer(kind=i8_kind))
554 call allocate_array(buf_i8_kind, e)
555 call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
556 unlim_dim_level=unlim_dim_level, &
557 corner=c, edge_lengths=e, broadcast=.false.)
561 if (buffer_includes_halos)
then
563 c(xdim_index) = isc - isd + 1
564 c(ydim_index) = jsc - jsd + 1
566 call put_array_section(buf_i8_kind, vdata, c, e)
569 call mpp_send(buf_i8_kind,
size(buf_i8_kind), fileobj%pelist(i))
572 deallocate(buf_i8_kind)
573 type is (real(kind=r4_kind))
575 call allocate_array(buf_r4_kind, e)
576 call netcdf_read_data(fileobj, variable_name, buf_r4_kind, &
577 unlim_dim_level=unlim_dim_level, &
578 corner=c, edge_lengths=e, broadcast=.false.)
582 if (buffer_includes_halos)
then
584 c(xdim_index) = isc - isd + 1
585 c(ydim_index) = jsc - jsd + 1
587 call put_array_section(buf_r4_kind, vdata, c, e)
590 call mpp_send(buf_r4_kind,
size(buf_r4_kind), fileobj%pelist(i))
593 deallocate(buf_r4_kind)
594 type is (real(kind=r8_kind))
596 call allocate_array(buf_r8_kind, e)
597 call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
598 unlim_dim_level=unlim_dim_level, &
599 corner=c, edge_lengths=e, broadcast=.false.)
603 if (buffer_includes_halos)
then
605 c(xdim_index) = isc - isd + 1
606 c(ydim_index) = jsc - jsd + 1
608 call put_array_section(buf_r8_kind, vdata, c, e)
611 call mpp_send(buf_r8_kind,
size(buf_r8_kind), fileobj%pelist(i))
614 deallocate(buf_r8_kind)
616 call error(
"unsupported variable type: domain_read_4d: file: "//trim(fileobj%path)//
" variable:"// &
617 & trim(variable_name))
621 deallocate(pe_icsize)
623 deallocate(pe_jcsize)
626 if (buffer_includes_halos)
then
627 c(xdim_index) = isc - isd + 1
628 c(ydim_index) = jsc - jsd + 1
630 e(xdim_index) = xc_size
631 e(ydim_index) = yc_size
633 type is (
integer(kind=i4_kind))
634 call allocate_array(buf_i4_kind, e)
635 call mpp_recv(buf_i4_kind,
size(buf_i4_kind), fileobj%io_root, block=.true.)
636 call put_array_section(buf_i4_kind, vdata, c, e)
637 deallocate(buf_i4_kind)
638 type is (
integer(kind=i8_kind))
639 call allocate_array(buf_i8_kind, e)
640 call mpp_recv(buf_i8_kind,
size(buf_i8_kind), fileobj%io_root, block=.true.)
641 call put_array_section(buf_i8_kind, vdata, c, e)
642 deallocate(buf_i8_kind)
643 type is (real(kind=r4_kind))
644 call allocate_array(buf_r4_kind, e)
645 call mpp_recv(buf_r4_kind,
size(buf_r4_kind), fileobj%io_root, block=.true.)
646 call put_array_section(buf_r4_kind, vdata, c, e)
647 deallocate(buf_r4_kind)
648 type is (real(kind=r8_kind))
649 call allocate_array(buf_r8_kind, e)
650 call mpp_recv(buf_r8_kind,
size(buf_r8_kind), fileobj%io_root, block=.true.)
651 call put_array_section(buf_r8_kind, vdata, c, e)
652 deallocate(buf_r8_kind)
654 call error(
"unsupported variable type: domain_read_4d: file: "//trim(fileobj%path)//
" variable:"// &
655 & trim(variable_name))
667 corner, edge_lengths)
669 type(fmsnetcdfdomainfile_t),
intent(in) :: fileobj
670 character(len=*),
intent(in) :: variable_name
671 class(*),
dimension(:,:,:,:,:),
intent(inout) :: vdata
674 integer,
intent(in),
optional :: unlim_dim_level
676 integer,
dimension(5),
intent(in),
optional :: corner
680 integer,
dimension(5),
intent(in),
optional :: edge_lengths
685 integer :: xdim_index
686 integer :: ydim_index
687 type(domain2d),
pointer :: io_domain
697 integer,
dimension(:),
allocatable :: pe_isc
698 integer,
dimension(:),
allocatable :: pe_icsize
699 integer,
dimension(:),
allocatable :: pe_jsc
700 integer,
dimension(:),
allocatable :: pe_jcsize
701 integer,
dimension(5) :: c
702 integer,
dimension(5) :: e
703 integer(kind=i4_kind),
dimension(:,:,:,:,:),
allocatable :: buf_i4_kind
704 integer(kind=i8_kind),
dimension(:,:,:,:,:),
allocatable :: buf_i8_kind
705 real(kind=r4_kind),
dimension(:,:,:,:,:),
allocatable :: buf_r4_kind
706 real(kind=r8_kind),
dimension(:,:,:,:,:),
allocatable :: buf_r8_kind
707 logical :: buffer_includes_halos
711 if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
712 xdim_index, ydim_index, xpos, ypos))
then
713 call netcdf_read_data(fileobj, variable_name, vdata, &
714 unlim_dim_level=unlim_dim_level, corner=corner, &
715 edge_lengths=edge_lengths, broadcast=.true.)
719 call domain_offsets(
size(vdata, xdim_index),
size(vdata, ydim_index), fileobj%domain, &
720 xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
721 msg=
"file:"//trim(fileobj%path)//
" and variable:"//trim(variable_name))
724 if (
present(edge_lengths)) e = edge_lengths
727 if (fileobj%is_root)
then
728 allocate(pe_isc(
size(fileobj%pelist)))
729 allocate(pe_icsize(
size(fileobj%pelist)))
730 allocate(pe_jsc(
size(fileobj%pelist)))
731 allocate(pe_jcsize(
size(fileobj%pelist)))
732 call mpp_get_compute_domains(io_domain, xbegin=pe_isc, xsize=pe_icsize, position=xpos)
733 call mpp_get_compute_domains(io_domain, ybegin=pe_jsc, ysize=pe_jcsize, position=ypos)
734 call mpp_get_global_domain(io_domain, xbegin=xgmin, position=xpos)
735 call mpp_get_global_domain(io_domain, ybegin=ygmin, position=ypos)
736 do i = 1,
size(fileobj%pelist)
738 if (
present(corner)) c = corner
739 c(xdim_index) = pe_isc(i)
740 c(ydim_index) = pe_jsc(i)
741 if (fileobj%adjust_indices)
then
742 c(xdim_index) = c(xdim_index) - xgmin + 1
743 c(ydim_index) = c(ydim_index) - ygmin + 1
745 e(xdim_index) = pe_icsize(i)
746 e(ydim_index) = pe_jcsize(i)
748 type is (
integer(kind=i4_kind))
750 call allocate_array(buf_i4_kind, e)
751 call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
752 unlim_dim_level=unlim_dim_level, &
753 corner=c, edge_lengths=e, broadcast=.false.)
758 if (buffer_includes_halos)
then
760 c(xdim_index) = isc - isd + 1
761 c(ydim_index) = jsc - jsd + 1
763 call put_array_section(buf_i4_kind, vdata, c, e)
766 call mpp_send(buf_i4_kind,
size(buf_i4_kind), fileobj%pelist(i))
769 deallocate(buf_i4_kind)
770 type is (
integer(kind=i8_kind))
772 call allocate_array(buf_i8_kind, e)
773 call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
774 unlim_dim_level=unlim_dim_level, &
775 corner=c, edge_lengths=e, broadcast=.false.)
779 if (buffer_includes_halos)
then
781 c(xdim_index) = isc - isd + 1
782 c(ydim_index) = jsc - jsd + 1
784 call put_array_section(buf_i8_kind, vdata, c, e)
787 call mpp_send(buf_i8_kind,
size(buf_i8_kind), fileobj%pelist(i))
790 deallocate(buf_i8_kind)
791 type is (real(kind=r4_kind))
793 call allocate_array(buf_r4_kind, e)
794 call netcdf_read_data(fileobj, variable_name, buf_r4_kind, &
795 unlim_dim_level=unlim_dim_level, &
796 corner=c, edge_lengths=e, broadcast=.false.)
800 if (buffer_includes_halos)
then
802 c(xdim_index) = isc - isd + 1
803 c(ydim_index) = jsc - jsd + 1
805 call put_array_section(buf_r4_kind, vdata, c, e)
808 call mpp_send(buf_r4_kind,
size(buf_r4_kind), fileobj%pelist(i))
811 deallocate(buf_r4_kind)
812 type is (real(kind=r8_kind))
814 call allocate_array(buf_r8_kind, e)
815 call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
816 unlim_dim_level=unlim_dim_level, &
817 corner=c, edge_lengths=e, broadcast=.false.)
821 if (buffer_includes_halos)
then
823 c(xdim_index) = isc - isd + 1
824 c(ydim_index) = jsc - jsd + 1
826 call put_array_section(buf_r8_kind, vdata, c, e)
829 call mpp_send(buf_r8_kind,
size(buf_r8_kind), fileobj%pelist(i))
832 deallocate(buf_r8_kind)
834 call error(
"unsupported variable type: domain_read_5d: file: "//trim(fileobj%path)//
" variable:"// &
835 & trim(variable_name))
839 deallocate(pe_icsize)
841 deallocate(pe_jcsize)
844 if (buffer_includes_halos)
then
845 c(xdim_index) = isc - isd + 1
846 c(ydim_index) = jsc - jsd + 1
848 e(xdim_index) = xc_size
849 e(ydim_index) = yc_size
851 type is (
integer(kind=i4_kind))
852 call allocate_array(buf_i4_kind, e)
853 call mpp_recv(buf_i4_kind,
size(buf_i4_kind), fileobj%io_root, block=.true.)
854 call put_array_section(buf_i4_kind, vdata, c, e)
855 deallocate(buf_i4_kind)
856 type is (
integer(kind=i8_kind))
857 call allocate_array(buf_i8_kind, e)
858 call mpp_recv(buf_i8_kind,
size(buf_i8_kind), fileobj%io_root, block=.true.)
859 call put_array_section(buf_i8_kind, vdata, c, e)
860 deallocate(buf_i8_kind)
861 type is (real(kind=r4_kind))
862 call allocate_array(buf_r4_kind, e)
863 call mpp_recv(buf_r4_kind,
size(buf_r4_kind), fileobj%io_root, block=.true.)
864 call put_array_section(buf_r4_kind, vdata, c, e)
865 deallocate(buf_r4_kind)
866 type is (real(kind=r8_kind))
867 call allocate_array(buf_r8_kind, e)
868 call mpp_recv(buf_r8_kind,
size(buf_r8_kind), fileobj%io_root, block=.true.)
869 call put_array_section(buf_r8_kind, vdata, c, e)
870 deallocate(buf_r8_kind)
872 call error(
"unsupported variable type: domain_read_5d: file: "//trim(fileobj%path)//
" variable:"// &
873 & trim(variable_name))
subroutine domain_read_0d(fileobj, variable_name, vdata, unlim_dim_level, corner)
I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and sca...
subroutine domain_read_3d(fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and sca...
subroutine domain_read_1d(fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and sca...
subroutine domain_read_2d(fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and sca...
subroutine domain_read_4d(fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and sca...
subroutine domain_read_5d(fileobj, variable_name, vdata, unlim_dim_level, corner, edge_lengths)
I/O domain root reads in a domain decomposed variable at a specific unlimited dimension level and sca...
type(domain2d) function, pointer mpp_get_io_domain(domain)
Set user stack size.
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...