30 subroutine domain_read_0d(fileobj, variable_name, vdata, unlim_dim_level, corner)
32 type(fmsnetcdfdomainfile_t),
intent(in) :: fileobj
33 character(len=*),
intent(in) :: variable_name
34 class(*),
intent(inout) :: vdata
37 integer,
intent(in),
optional :: unlim_dim_level
39 integer,
intent(in),
optional :: corner
44 call netcdf_read_data(fileobj, variable_name, vdata, &
45 unlim_dim_level=unlim_dim_level, corner=corner, &
59 type(fmsnetcdfdomainfile_t),
intent(in) :: fileobj
60 character(len=*),
intent(in) :: variable_name
61 class(*),
dimension(:),
intent(inout) :: vdata
64 integer,
intent(in),
optional :: unlim_dim_level
66 integer,
dimension(1),
intent(in),
optional :: corner
70 integer,
dimension(1),
intent(in),
optional :: edge_lengths
75 call netcdf_read_data(fileobj, variable_name, vdata, &
76 unlim_dim_level=unlim_dim_level, corner=corner, &
77 edge_lengths=edge_lengths, broadcast=.true.)
89 type(fmsnetcdfdomainfile_t),
intent(in) :: fileobj
90 character(len=*),
intent(in) :: variable_name
91 class(*),
contiguous,
target,
intent(inout) :: vdata(:,:)
94 integer,
intent(in),
optional :: unlim_dim_level
96 integer,
dimension(2),
intent(in),
optional :: corner
100 integer,
dimension(2),
intent(in),
optional :: edge_lengths
105 integer :: xdim_index
106 integer :: ydim_index
118 logical :: buffer_includes_halos
123 type(domain2d),
pointer :: io_domain
126 integer(kind=i4_kind),
dimension(:,:),
allocatable :: buf_i4_kind_pe
127 integer(kind=i8_kind),
dimension(:,:),
allocatable :: buf_i8_kind_pe
128 real(kind=r4_kind),
dimension(:,:),
allocatable :: buf_r4_kind_pe
129 real(kind=r8_kind),
dimension(:,:),
allocatable :: buf_r8_kind_pe
130 integer(kind=i4_kind),
dimension(:,:),
allocatable :: buf_i4_kind
131 integer(kind=i8_kind),
dimension(:,:),
allocatable :: buf_i8_kind
132 real(kind=r4_kind),
dimension(:,:),
allocatable :: buf_r4_kind
133 real(kind=r8_kind),
dimension(:,:),
allocatable :: buf_r8_kind
134 class(*),
dimension(:,:,:,:),
pointer :: vdata_dummy
136 if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
137 xdim_index, ydim_index, xpos, ypos))
then
138 call netcdf_read_data(fileobj, variable_name, vdata, &
139 unlim_dim_level=unlim_dim_level, corner=corner, &
140 edge_lengths=edge_lengths, broadcast=.true.)
144 if (xdim_index .ne. 1 .or. ydim_index .ne. 2)
then
148 vdata_dummy(1:
size(vdata,1),1:
size(vdata,2), 1:1, 1:1) => vdata(:,:)
149 call domain_read_4d(fileobj, variable_name, vdata_dummy, unlim_dim_level)
153 call domain_offsets(
size(vdata, xdim_index),
size(vdata, ydim_index), fileobj%domain, &
154 xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
155 msg=
"file:"//trim(fileobj%path)//
" and variable:"//trim(variable_name))
159 call mpp_get_global_domain(io_domain, xbegin=xgbegin, xsize=xgsize, position=xpos)
160 call mpp_get_global_domain(io_domain, ybegin=ygbegin, ysize=ygsize, position=ypos)
163 if (fileobj%is_root)
then
165 if (fileobj%adjust_indices)
then
172 c(xdim_index) = xgbegin
173 c(ydim_index) = ygbegin
176 e(xdim_index) = xgsize
177 e(ydim_index) = ygsize
181 type is (
integer(kind=i4_kind))
182 call allocate_array(buf_i4_kind, e)
183 call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
184 unlim_dim_level=unlim_dim_level, &
185 corner=c, edge_lengths=e, broadcast=.false.)
186 type is (
integer(kind=i8_kind))
187 call allocate_array(buf_i8_kind, e)
188 call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
189 unlim_dim_level=unlim_dim_level, &
190 corner=c, edge_lengths=e, broadcast=.false.)
191 type is (real(kind=r4_kind))
192 call allocate_array(buf_r4_kind, e)
193 call netcdf_read_data(fileobj, variable_name, buf_r4_kind, &
194 unlim_dim_level=unlim_dim_level, &
195 corner=c, edge_lengths=e, broadcast=.false.)
196 type is (real(kind=r8_kind))
197 call allocate_array(buf_r8_kind, e)
198 call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
199 unlim_dim_level=unlim_dim_level, &
200 corner=c, edge_lengths=e, broadcast=.false.)
202 call error(
"unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//
" variable:"// &
203 & trim(variable_name))
211 if (buffer_includes_halos)
then
213 c(xdim_index) = isc - isd + 1
214 c(ydim_index) = jsc - jsd + 1
220 e(xdim_index) = xc_size
221 e(ydim_index) = yc_size
224 type is (
integer(kind=i4_kind))
225 call allocate_array(buf_i4_kind_pe, e)
226 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, fileobj%pelist, &
227 buf_i4_kind_pe, buf_i4_kind, fileobj%is_root)
228 call put_array_section(buf_i4_kind_pe, vdata, c, e)
229 deallocate(buf_i4_kind_pe)
230 type is (
integer(kind=i8_kind))
231 call allocate_array(buf_i8_kind_pe, e)
232 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, fileobj%pelist, &
233 buf_i8_kind_pe, buf_i8_kind, fileobj%is_root)
234 call put_array_section(buf_i8_kind_pe, vdata, c, e)
235 deallocate(buf_i8_kind_pe)
236 type is (real(kind=r4_kind))
237 call allocate_array(buf_r4_kind_pe, e)
238 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, fileobj%pelist, &
239 buf_r4_kind_pe, buf_r4_kind, fileobj%is_root)
240 call put_array_section(buf_r4_kind_pe, vdata, c, e)
241 deallocate(buf_r4_kind_pe)
242 type is (real(kind=r8_kind))
243 call allocate_array(buf_r8_kind_pe, e)
244 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, fileobj%pelist, &
245 buf_r8_kind_pe, buf_r8_kind, fileobj%is_root)
246 call put_array_section(buf_r8_kind_pe, vdata, c, e)
247 deallocate(buf_r8_kind_pe)
249 call error(
"unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//
" variable:"// &
250 & trim(variable_name))
253 if (fileobj%is_root)
then
254 if (
allocated(buf_i4_kind))
deallocate(buf_i4_kind)
255 if (
allocated(buf_i8_kind))
deallocate(buf_i8_kind)
256 if (
allocated(buf_r4_kind))
deallocate(buf_r4_kind)
257 if (
allocated(buf_r8_kind))
deallocate(buf_r8_kind)
268 corner, edge_lengths)
269 type(fmsnetcdfdomainfile_t),
intent(in) :: fileobj
270 character(len=*),
intent(in) :: variable_name
271 class(*),
contiguous,
target,
intent(inout) :: vdata(:,:,:)
274 integer,
intent(in),
optional :: unlim_dim_level
276 integer,
dimension(3),
intent(in),
optional :: corner
280 integer,
dimension(3),
intent(in),
optional :: edge_lengths
285 integer :: xdim_index
286 integer :: ydim_index
298 logical :: buffer_includes_halos
303 type(domain2d),
pointer :: io_domain
306 integer(kind=i4_kind),
dimension(:,:,:),
allocatable :: buf_i4_kind_pe
307 integer(kind=i8_kind),
dimension(:,:,:),
allocatable :: buf_i8_kind_pe
308 real(kind=r4_kind),
dimension(:,:,:),
allocatable :: buf_r4_kind_pe
309 real(kind=r8_kind),
dimension(:,:,:),
allocatable :: buf_r8_kind_pe
310 integer(kind=i4_kind),
dimension(:,:,:),
allocatable :: buf_i4_kind
311 integer(kind=i8_kind),
dimension(:,:,:),
allocatable :: buf_i8_kind
312 real(kind=r4_kind),
dimension(:,:,:),
allocatable :: buf_r4_kind
313 real(kind=r8_kind),
dimension(:,:,:),
allocatable :: buf_r8_kind
314 class(*),
dimension(:,:,:,:),
pointer :: vdata_dummy
316 if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
317 xdim_index, ydim_index, xpos, ypos))
then
318 call netcdf_read_data(fileobj, variable_name, vdata, &
319 unlim_dim_level=unlim_dim_level, corner=corner, &
320 edge_lengths=edge_lengths, broadcast=.true.)
324 if (xdim_index .ne. 1 .or. ydim_index .ne. 2)
then
328 vdata_dummy(1:
size(vdata,1),1:
size(vdata,2), 1:
size(vdata,3), 1:1) => vdata(:,:,:)
329 call domain_read_4d(fileobj, variable_name, vdata_dummy, unlim_dim_level)
333 call domain_offsets(
size(vdata, xdim_index),
size(vdata, ydim_index), fileobj%domain, &
334 xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
335 msg=
"file:"//trim(fileobj%path)//
" and variable:"//trim(variable_name))
339 call mpp_get_global_domain(io_domain, xbegin=xgbegin, xsize=xgsize, position=xpos)
340 call mpp_get_global_domain(io_domain, ybegin=ygbegin, ysize=ygsize, position=ypos)
343 if (fileobj%is_root)
then
345 if (fileobj%adjust_indices)
then
352 c(xdim_index) = xgbegin
353 c(ydim_index) = ygbegin
356 e(xdim_index) = xgsize
357 e(ydim_index) = ygsize
361 type is (
integer(kind=i4_kind))
362 call allocate_array(buf_i4_kind, e)
363 call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
364 unlim_dim_level=unlim_dim_level, &
365 corner=c, edge_lengths=e, broadcast=.false.)
366 type is (
integer(kind=i8_kind))
367 call allocate_array(buf_i8_kind, e)
368 call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
369 unlim_dim_level=unlim_dim_level, &
370 corner=c, edge_lengths=e, broadcast=.false.)
371 type is (real(kind=r4_kind))
372 call allocate_array(buf_r4_kind, e)
373 call netcdf_read_data(fileobj, variable_name, buf_r4_kind, &
374 unlim_dim_level=unlim_dim_level, &
375 corner=c, edge_lengths=e, broadcast=.false.)
376 type is (real(kind=r8_kind))
377 call allocate_array(buf_r8_kind, e)
378 call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
379 unlim_dim_level=unlim_dim_level, &
380 corner=c, edge_lengths=e, broadcast=.false.)
382 call error(
"unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//
" variable:"// &
383 & trim(variable_name))
391 if (buffer_includes_halos)
then
393 c(xdim_index) = isc - isd + 1
394 c(ydim_index) = jsc - jsd + 1
400 e(xdim_index) = xc_size
401 e(ydim_index) = yc_size
404 type is (
integer(kind=i4_kind))
405 call allocate_array(buf_i4_kind_pe, e)
406 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
407 buf_i4_kind_pe, buf_i4_kind, fileobj%is_root)
408 call put_array_section(buf_i4_kind_pe, vdata, c, e)
409 deallocate(buf_i4_kind_pe)
410 type is (
integer(kind=i8_kind))
411 call allocate_array(buf_i8_kind_pe, e)
412 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
413 buf_i8_kind_pe, buf_i8_kind, fileobj%is_root)
414 call put_array_section(buf_i8_kind_pe, vdata, c, e)
415 deallocate(buf_i8_kind_pe)
416 type is (real(kind=r4_kind))
417 call allocate_array(buf_r4_kind_pe, e)
418 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
419 buf_r4_kind_pe, buf_r4_kind, fileobj%is_root)
420 call put_array_section(buf_r4_kind_pe, vdata, c, e)
421 deallocate(buf_r4_kind_pe)
422 type is (real(kind=r8_kind))
423 call allocate_array(buf_r8_kind_pe, e)
424 call mpp_scatter(isc-xgbegin+1, isc+xc_size-xgbegin, jsc-ygbegin+1, jsc+yc_size-ygbegin, e(3), fileobj%pelist, &
425 buf_r8_kind_pe, buf_r8_kind, fileobj%is_root)
426 call put_array_section(buf_r8_kind_pe, vdata, c, e)
427 deallocate(buf_r8_kind_pe)
429 call error(
"unsupported variable type: domain_read_2d: file: "//trim(fileobj%path)//
" variable:"// &
430 & trim(variable_name))
433 if (fileobj%is_root)
then
434 if (
allocated(buf_i4_kind))
deallocate(buf_i4_kind)
435 if (
allocated(buf_i8_kind))
deallocate(buf_i8_kind)
436 if (
allocated(buf_r4_kind))
deallocate(buf_r4_kind)
437 if (
allocated(buf_r8_kind))
deallocate(buf_r8_kind)
449 corner, edge_lengths)
451 type(fmsnetcdfdomainfile_t),
intent(in) :: fileobj
452 character(len=*),
intent(in) :: variable_name
453 class(*),
dimension(:,:,:,:),
intent(inout) :: vdata
456 integer,
intent(in),
optional :: unlim_dim_level
458 integer,
dimension(4),
intent(in),
optional :: corner
462 integer,
dimension(4),
intent(in),
optional :: edge_lengths
467 integer :: xdim_index
468 integer :: ydim_index
469 type(domain2d),
pointer :: io_domain
479 integer,
dimension(:),
allocatable :: pe_isc
480 integer,
dimension(:),
allocatable :: pe_icsize
481 integer,
dimension(:),
allocatable :: pe_jsc
482 integer,
dimension(:),
allocatable :: pe_jcsize
483 integer,
dimension(4) :: c
484 integer,
dimension(4) :: e
485 integer(kind=i4_kind),
dimension(:,:,:,:),
allocatable :: buf_i4_kind
486 integer(kind=i8_kind),
dimension(:,:,:,:),
allocatable :: buf_i8_kind
487 real(kind=r4_kind),
dimension(:,:,:,:),
allocatable :: buf_r4_kind
488 real(kind=r8_kind),
dimension(:,:,:,:),
allocatable :: buf_r8_kind
489 logical :: buffer_includes_halos
493 if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
494 xdim_index, ydim_index, xpos, ypos))
then
495 call netcdf_read_data(fileobj, variable_name, vdata, &
496 unlim_dim_level=unlim_dim_level, corner=corner, &
497 edge_lengths=edge_lengths, broadcast=.true.)
501 call domain_offsets(
size(vdata, xdim_index),
size(vdata, ydim_index), fileobj%domain, &
502 xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
503 msg=
"file:"//trim(fileobj%path)//
" and variable:"//trim(variable_name))
508 if (fileobj%is_root)
then
509 allocate(pe_isc(
size(fileobj%pelist)))
510 allocate(pe_icsize(
size(fileobj%pelist)))
511 allocate(pe_jsc(
size(fileobj%pelist)))
512 allocate(pe_jcsize(
size(fileobj%pelist)))
513 call mpp_get_compute_domains(io_domain, xbegin=pe_isc, xsize=pe_icsize, position=xpos)
514 call mpp_get_compute_domains(io_domain, ybegin=pe_jsc, ysize=pe_jcsize, position=ypos)
515 call mpp_get_global_domain(io_domain, xbegin=xgmin, position=xpos)
516 call mpp_get_global_domain(io_domain, ybegin=ygmin, position=ypos)
517 do i = 1,
size(fileobj%pelist)
518 c(xdim_index) = pe_isc(i)
519 c(ydim_index) = pe_jsc(i)
520 if (fileobj%adjust_indices)
then
521 c(xdim_index) = c(xdim_index) - xgmin + 1
522 c(ydim_index) = c(ydim_index) - ygmin + 1
524 e(xdim_index) = pe_icsize(i)
525 e(ydim_index) = pe_jcsize(i)
527 type is (
integer(kind=i4_kind))
529 call allocate_array(buf_i4_kind, e)
530 call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
531 unlim_dim_level=unlim_dim_level, &
532 corner=c, edge_lengths=e, broadcast=.false.)
535 if (buffer_includes_halos)
then
537 c(xdim_index) = isc - isd + 1
538 c(ydim_index) = jsc - jsd + 1
543 call put_array_section(buf_i4_kind, vdata, c, e)
546 call mpp_send(buf_i4_kind,
size(buf_i4_kind), fileobj%pelist(i))
549 deallocate(buf_i4_kind)
550 type is (
integer(kind=i8_kind))
552 call allocate_array(buf_i8_kind, e)
553 call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
554 unlim_dim_level=unlim_dim_level, &
555 corner=c, edge_lengths=e, broadcast=.false.)
558 if (buffer_includes_halos)
then
560 c(xdim_index) = isc - isd + 1
561 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.)
581 if (buffer_includes_halos)
then
583 c(xdim_index) = isc - isd + 1
584 c(ydim_index) = jsc - jsd + 1
589 call put_array_section(buf_r4_kind, vdata, c, e)
592 call mpp_send(buf_r4_kind,
size(buf_r4_kind), fileobj%pelist(i))
595 deallocate(buf_r4_kind)
596 type is (real(kind=r8_kind))
598 call allocate_array(buf_r8_kind, e)
599 call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
600 unlim_dim_level=unlim_dim_level, &
601 corner=c, edge_lengths=e, broadcast=.false.)
604 if (buffer_includes_halos)
then
606 c(xdim_index) = isc - isd + 1
607 c(ydim_index) = jsc - jsd + 1
612 call put_array_section(buf_r8_kind, vdata, c, e)
615 call mpp_send(buf_r8_kind,
size(buf_r8_kind), fileobj%pelist(i))
618 deallocate(buf_r8_kind)
620 call error(
"unsupported variable type: domain_read_4d: file: "//trim(fileobj%path)//
" variable:"// &
621 & trim(variable_name))
625 deallocate(pe_icsize)
627 deallocate(pe_jcsize)
629 if (buffer_includes_halos)
then
630 c(xdim_index) = isc - isd + 1
631 c(ydim_index) = jsc - jsd + 1
633 e(xdim_index) = xc_size
634 e(ydim_index) = yc_size
636 type is (
integer(kind=i4_kind))
637 call allocate_array(buf_i4_kind, e)
638 call mpp_recv(buf_i4_kind,
size(buf_i4_kind), fileobj%io_root, block=.true.)
639 call put_array_section(buf_i4_kind, vdata, c, e)
640 deallocate(buf_i4_kind)
641 type is (
integer(kind=i8_kind))
642 call allocate_array(buf_i8_kind, e)
643 call mpp_recv(buf_i8_kind,
size(buf_i8_kind), fileobj%io_root, block=.true.)
644 call put_array_section(buf_i8_kind, vdata, c, e)
645 deallocate(buf_i8_kind)
646 type is (real(kind=r4_kind))
647 call allocate_array(buf_r4_kind, e)
648 call mpp_recv(buf_r4_kind,
size(buf_r4_kind), fileobj%io_root, block=.true.)
649 call put_array_section(buf_r4_kind, vdata, c, e)
650 deallocate(buf_r4_kind)
651 type is (real(kind=r8_kind))
652 call allocate_array(buf_r8_kind, e)
653 call mpp_recv(buf_r8_kind,
size(buf_r8_kind), fileobj%io_root, block=.true.)
654 call put_array_section(buf_r8_kind, vdata, c, e)
655 deallocate(buf_r8_kind)
657 call error(
"unsupported variable type: domain_read_4d: file: "//trim(fileobj%path)//
" variable:"// &
658 & trim(variable_name))
670 corner, edge_lengths)
672 type(fmsnetcdfdomainfile_t),
intent(in) :: fileobj
673 character(len=*),
intent(in) :: variable_name
674 class(*),
dimension(:,:,:,:,:),
intent(inout) :: vdata
677 integer,
intent(in),
optional :: unlim_dim_level
679 integer,
dimension(5),
intent(in),
optional :: corner
683 integer,
dimension(5),
intent(in),
optional :: edge_lengths
688 integer :: xdim_index
689 integer :: ydim_index
690 type(domain2d),
pointer :: io_domain
700 integer,
dimension(:),
allocatable :: pe_isc
701 integer,
dimension(:),
allocatable :: pe_icsize
702 integer,
dimension(:),
allocatable :: pe_jsc
703 integer,
dimension(:),
allocatable :: pe_jcsize
704 integer,
dimension(5) :: c
705 integer,
dimension(5) :: e
706 integer(kind=i4_kind),
dimension(:,:,:,:,:),
allocatable :: buf_i4_kind
707 integer(kind=i8_kind),
dimension(:,:,:,:,:),
allocatable :: buf_i8_kind
708 real(kind=r4_kind),
dimension(:,:,:,:,:),
allocatable :: buf_r4_kind
709 real(kind=r8_kind),
dimension(:,:,:,:,:),
allocatable :: buf_r8_kind
710 logical :: buffer_includes_halos
714 if (.not. is_variable_domain_decomposed(fileobj, variable_name, .true., &
715 xdim_index, ydim_index, xpos, ypos))
then
716 call netcdf_read_data(fileobj, variable_name, vdata, &
717 unlim_dim_level=unlim_dim_level, corner=corner, &
718 edge_lengths=edge_lengths, broadcast=.true.)
722 call domain_offsets(
size(vdata, xdim_index),
size(vdata, ydim_index), fileobj%domain, &
723 xpos, ypos, isd, isc, xc_size, jsd, jsc, yc_size, buffer_includes_halos, &
724 msg=
"file:"//trim(fileobj%path)//
" and variable:"//trim(variable_name))
729 if (fileobj%is_root)
then
730 allocate(pe_isc(
size(fileobj%pelist)))
731 allocate(pe_icsize(
size(fileobj%pelist)))
732 allocate(pe_jsc(
size(fileobj%pelist)))
733 allocate(pe_jcsize(
size(fileobj%pelist)))
734 call mpp_get_compute_domains(io_domain, xbegin=pe_isc, xsize=pe_icsize, position=xpos)
735 call mpp_get_compute_domains(io_domain, ybegin=pe_jsc, ysize=pe_jcsize, position=ypos)
736 call mpp_get_global_domain(io_domain, xbegin=xgmin, position=xpos)
737 call mpp_get_global_domain(io_domain, ybegin=ygmin, position=ypos)
738 do i = 1,
size(fileobj%pelist)
740 c(xdim_index) = pe_isc(i)
741 c(ydim_index) = pe_jsc(i)
742 if (fileobj%adjust_indices)
then
743 c(xdim_index) = c(xdim_index) - xgmin + 1
744 c(ydim_index) = c(ydim_index) - ygmin + 1
746 e(xdim_index) = pe_icsize(i)
747 e(ydim_index) = pe_jcsize(i)
749 type is (
integer(kind=i4_kind))
751 call allocate_array(buf_i4_kind, e)
752 call netcdf_read_data(fileobj, variable_name, buf_i4_kind, &
753 unlim_dim_level=unlim_dim_level, &
754 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
766 call put_array_section(buf_i4_kind, vdata, c, e)
769 call mpp_send(buf_i4_kind,
size(buf_i4_kind), fileobj%pelist(i))
772 deallocate(buf_i4_kind)
773 type is (
integer(kind=i8_kind))
775 call allocate_array(buf_i8_kind, e)
776 call netcdf_read_data(fileobj, variable_name, buf_i8_kind, &
777 unlim_dim_level=unlim_dim_level, &
778 corner=c, edge_lengths=e, broadcast=.false.)
781 if (buffer_includes_halos)
then
783 c(xdim_index) = isc - isd + 1
784 c(ydim_index) = jsc - jsd + 1
789 call put_array_section(buf_i8_kind, vdata, c, e)
792 call mpp_send(buf_i8_kind,
size(buf_i8_kind), fileobj%pelist(i))
795 deallocate(buf_i8_kind)
796 type is (real(kind=r4_kind))
798 call allocate_array(buf_r4_kind, e)
799 call netcdf_read_data(fileobj, variable_name, buf_r4_kind, &
800 unlim_dim_level=unlim_dim_level, &
801 corner=c, edge_lengths=e, broadcast=.false.)
804 if (buffer_includes_halos)
then
806 c(xdim_index) = isc - isd + 1
807 c(ydim_index) = jsc - jsd + 1
812 call put_array_section(buf_r4_kind, vdata, c, e)
815 call mpp_send(buf_r4_kind,
size(buf_r4_kind), fileobj%pelist(i))
818 deallocate(buf_r4_kind)
819 type is (real(kind=r8_kind))
821 call allocate_array(buf_r8_kind, e)
822 call netcdf_read_data(fileobj, variable_name, buf_r8_kind, &
823 unlim_dim_level=unlim_dim_level, &
824 corner=c, edge_lengths=e, broadcast=.false.)
827 if (buffer_includes_halos)
then
829 c(xdim_index) = isc - isd + 1
830 c(ydim_index) = jsc - jsd + 1
835 call put_array_section(buf_r8_kind, vdata, c, e)
838 call mpp_send(buf_r8_kind,
size(buf_r8_kind), fileobj%pelist(i))
841 deallocate(buf_r8_kind)
843 call error(
"unsupported variable type: domain_read_5d: file: "//trim(fileobj%path)//
" variable:"// &
844 & trim(variable_name))
848 deallocate(pe_icsize)
850 deallocate(pe_jcsize)
852 if (buffer_includes_halos)
then
853 c(xdim_index) = isc - isd + 1
854 c(ydim_index) = jsc - jsd + 1
856 e(xdim_index) = xc_size
857 e(ydim_index) = yc_size
859 type is (
integer(kind=i4_kind))
860 call allocate_array(buf_i4_kind, e)
861 call mpp_recv(buf_i4_kind,
size(buf_i4_kind), fileobj%io_root, block=.true.)
862 call put_array_section(buf_i4_kind, vdata, c, e)
863 deallocate(buf_i4_kind)
864 type is (
integer(kind=i8_kind))
865 call allocate_array(buf_i8_kind, e)
866 call mpp_recv(buf_i8_kind,
size(buf_i8_kind), fileobj%io_root, block=.true.)
867 call put_array_section(buf_i8_kind, vdata, c, e)
868 deallocate(buf_i8_kind)
869 type is (real(kind=r4_kind))
870 call allocate_array(buf_r4_kind, e)
871 call mpp_recv(buf_r4_kind,
size(buf_r4_kind), fileobj%io_root, block=.true.)
872 call put_array_section(buf_r4_kind, vdata, c, e)
873 deallocate(buf_r4_kind)
874 type is (real(kind=r8_kind))
875 call allocate_array(buf_r8_kind, e)
876 call mpp_recv(buf_r8_kind,
size(buf_r8_kind), fileobj%io_root, block=.true.)
877 call put_array_section(buf_r8_kind, vdata, c, e)
878 deallocate(buf_r8_kind)
880 call error(
"unsupported variable type: domain_read_5d: file: "//trim(fileobj%path)//
" variable:"// &
881 & 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...