26 subroutine scatter_data_bc_2d(fileobj, varname, vdata, bc_info, unlim_dim_level, ignore_checksum)
27 class(fmsnetcdffile_t),
intent(inout) :: fileobj
28 character(len=*),
intent(in) :: varname
29 class(*),
dimension(:,:),
intent(inout) :: vdata
30 type(bc_information),
intent(inout) :: bc_info
31 integer,
intent(in),
optional :: unlim_dim_level
32 logical,
intent(in),
optional :: ignore_checksum
34 real(kind=r8_kind),
dimension(:,:),
allocatable :: global_buf_r8_kind
35 real(kind=r8_kind),
dimension(:,:),
allocatable :: local_buf_r8_kind
36 real(kind=r4_kind),
dimension(:,:),
allocatable :: global_buf_r4_kind
37 real(kind=r4_kind),
dimension(:,:),
allocatable :: local_buf_r4_kind
39 integer(kind=i8_kind) :: chksum_val
40 character(len=32) :: chksum
41 character(len=32) :: chksum_read
42 logical :: chksum_ignore = .false.
45 integer :: isc, iec, jsc, jec
46 integer :: i1, i2, j1, j2
47 integer :: i_add, j_add
52 if (
PRESENT(ignore_checksum)) chksum_ignore = ignore_checksum
55 isc = bc_info%indices(1)
56 iec = bc_info%indices(2)
57 jsc = bc_info%indices(3)
58 jec = bc_info%indices(4)
61 i1 = 1 + bc_info%x_halo
63 j1 = 1 + bc_info%y_halo
67 i_add = bc_info%ishift
68 j_add = bc_info%jshift
70 if (fileobj%is_root)
then
71 i_glob = bc_info%global_size(1)
72 j_glob = bc_info%global_size(2)
76 type is (real(kind=r4_kind))
77 if (fileobj%is_root)
then
79 allocate(global_buf_r4_kind(i_glob, j_glob))
81 call netcdf_read_data(fileobj, varname, &
83 unlim_dim_level=unlim_dim_level, &
86 if (.not.chksum_ignore .and. variable_att_exists(fileobj, varname,
"checksum", broadcast=.false.))
then
87 call get_variable_attribute(fileobj, varname,
"checksum", chksum_read, broadcast=.false.)
88 chksum_val = mpp_chksum(global_buf_r4_kind, (/
mpp_pe()/))
90 write(chksum,
"(Z16)") chksum_val
91 if (.not. string_compare(adjustl(chksum), adjustl(chksum_read)))
then
92 call error(
"scatter_data_bc_2d: "//trim(varname)//
" chksum_in:"//(chksum)//
" chksum_out:"// &
100 if (fileobj%is_root .and. .not. bc_info%data_on_file_root)
then
101 allocate(local_buf_r4_kind(1,1))
102 local_buf_r4_kind = 0.
103 isc = 1; i_add=0; iec=1; jsc=1; j_add=0; jec=1
104 i1=1; i2=1; j1=1; j2=1
106 allocate(local_buf_r4_kind(i1:i2,j1:j2))
110 call mpp_scatter(isc+i_add, iec+i_add, jsc+j_add, jec+j_add, bc_info%pelist, &
111 local_buf_r4_kind(i1:i2,j1:j2), global_buf_r4_kind, fileobj%is_root)
114 if (fileobj%is_root .and. .not. bc_info%data_on_file_root)
then
115 deallocate(local_buf_r4_kind)
120 vdata(i1:i2,j1:j2) = local_buf_r4_kind(i1:i2,j1:j2)
122 deallocate(local_buf_r4_kind)
123 if (fileobj%is_root)
deallocate(global_buf_r4_kind)
125 type is (real(kind=r8_kind))
126 if (fileobj%is_root)
then
128 allocate(global_buf_r8_kind(i_glob, j_glob))
130 call netcdf_read_data(fileobj, varname, &
131 global_buf_r8_kind, &
132 unlim_dim_level=unlim_dim_level, &
135 if (.not.chksum_ignore .and. variable_att_exists(fileobj, varname,
"checksum", broadcast=.false.))
then
136 call get_variable_attribute(fileobj, varname,
"checksum", chksum_read, broadcast=.false.)
137 chksum_val = mpp_chksum(global_buf_r8_kind, (/
mpp_pe()/))
139 write(chksum,
"(Z16)") chksum_val
140 if (.not. string_compare(adjustl(chksum), adjustl(chksum_read)))
then
141 call error(
"scatter_data_bc_2d: "//trim(varname)//
" chksum_in:"//chksum//
" chksum_out:"//chksum_read)
148 if (fileobj%is_root .and. .not. bc_info%data_on_file_root)
then
149 allocate(local_buf_r8_kind(1,1))
150 local_buf_r8_kind = 0.
151 isc = 1; i_add=0; iec=1; jsc=1; j_add=0; jec=1
152 i1=1; i2=1; j1=1; j2=1
154 allocate(local_buf_r8_kind(i1:i2,j1:j2))
158 call mpp_scatter(isc+i_add, iec+i_add, jsc+j_add, jec+j_add, bc_info%pelist, &
159 local_buf_r8_kind(i1:i2,j1:j2), global_buf_r8_kind, fileobj%is_root)
162 if (fileobj%is_root .and. .not. bc_info%data_on_file_root)
then
163 deallocate(local_buf_r8_kind)
168 vdata(i1:i2,j1:j2) = local_buf_r8_kind(i1:i2,j1:j2)
170 deallocate(local_buf_r8_kind)
171 if (fileobj%is_root)
deallocate(global_buf_r8_kind)
174 call error(
"scatter_data_bc_2d: unsupported type. Currently only r8_kind and r4_kinds are supported")
180 class(fmsnetcdffile_t),
intent(inout) :: fileobj
181 character(len=*),
intent(in) :: varname
182 class(*),
dimension(:,:,:),
intent(inout) :: vdata
183 type(bc_information),
intent(inout) :: bc_info
184 integer,
intent(in),
optional :: unlim_dim_level
185 logical,
intent(in),
optional :: ignore_checksum
187 real(kind=r8_kind),
dimension(:,:,:),
allocatable :: global_buf_r8_kind
188 real(kind=r8_kind),
dimension(:,:,:),
allocatable :: local_buf_r8_kind
189 real(kind=r4_kind),
dimension(:,:,:),
allocatable :: global_buf_r4_kind
190 real(kind=r4_kind),
dimension(:,:,:),
allocatable :: local_buf_r4_kind
192 integer(kind=i8_kind) :: chksum_val
193 character(len=32) :: chksum
194 character(len=32) :: chksum_read
195 logical :: chksum_ignore = .false.
198 integer :: isc, iec, jsc, jec
199 integer :: i1, i2, j1, j2
200 integer :: i_add, j_add
206 if (
PRESENT(ignore_checksum)) chksum_ignore = ignore_checksum
209 isc = bc_info%indices(1)
210 iec = bc_info%indices(2)
211 jsc = bc_info%indices(3)
212 jec = bc_info%indices(4)
215 i1 = 1 + bc_info%x_halo
217 j1 = 1 + bc_info%y_halo
221 i_add = bc_info%ishift
222 j_add = bc_info%jshift
224 if (fileobj%is_root)
then
225 i_glob = bc_info%global_size(1)
226 j_glob = bc_info%global_size(2)
229 k_glob = bc_info%global_size(3)
232 type is (real(kind=r4_kind))
233 if (fileobj%is_root)
then
235 allocate(global_buf_r4_kind(i_glob, j_glob, k_glob))
237 call netcdf_read_data(fileobj, varname, &
238 global_buf_r4_kind, &
239 unlim_dim_level=unlim_dim_level, &
242 if (.not.chksum_ignore .and. variable_att_exists(fileobj, varname,
"checksum", broadcast=.false.))
then
243 call get_variable_attribute(fileobj, varname,
"checksum", chksum_read, broadcast=.false.)
244 chksum_val = mpp_chksum(global_buf_r4_kind, (/
mpp_pe()/))
246 write(chksum,
"(Z16)") chksum_val
247 if (.not. string_compare(adjustl(chksum), adjustl(chksum_read)))
then
248 call error(
"scatter_data_bc_2d: "//trim(varname)//
" chksum_in:"//(chksum)//
" chksum_out:"// &
256 if (fileobj%is_root .and. .not. bc_info%data_on_file_root)
then
257 allocate(local_buf_r4_kind(1,1,k_glob))
258 local_buf_r4_kind = 0.
259 isc = 1; i_add=0; iec=1; jsc=1; j_add=0; jec=1
260 i1=1; i2=1; j1=1; j2=1
262 allocate(local_buf_r4_kind(i1:i2,j1:j2,k_glob))
266 call mpp_scatter(isc+i_add, iec+i_add, jsc+j_add, jec+j_add, k_glob, bc_info%pelist, &
267 local_buf_r4_kind(i1:i2,j1:j2,:), global_buf_r4_kind, fileobj%is_root)
270 if (fileobj%is_root .and. .not. bc_info%data_on_file_root)
then
271 deallocate(local_buf_r4_kind)
276 vdata(i1:i2,j1:j2,:) = local_buf_r4_kind(i1:i2,j1:j2,:)
278 deallocate(local_buf_r4_kind)
279 if (fileobj%is_root)
deallocate(global_buf_r4_kind)
281 type is (real(kind=r8_kind))
282 if (fileobj%is_root)
then
284 allocate(global_buf_r8_kind(i_glob, j_glob, k_glob))
286 call netcdf_read_data(fileobj, varname, &
287 global_buf_r8_kind, &
288 unlim_dim_level=unlim_dim_level, &
291 if (.not.chksum_ignore .and. variable_att_exists(fileobj, varname,
"checksum", broadcast=.false.))
then
292 call get_variable_attribute(fileobj, varname,
"checksum", chksum_read, broadcast=.false.)
293 chksum_val = mpp_chksum(global_buf_r8_kind, (/
mpp_pe()/))
295 write(chksum,
"(Z16)") chksum_val
296 if (.not. string_compare(adjustl(chksum), adjustl(chksum_read)))
then
297 call error(
"scatter_data_bc_2d: "//trim(varname)//
" chksum_in:"//chksum//
" chksum_out:"//chksum_read)
304 if (fileobj%is_root .and. .not. bc_info%data_on_file_root)
then
305 allocate(local_buf_r8_kind(1,1,k_glob))
306 local_buf_r8_kind = 0.
307 isc = 1; i_add=0; iec=1; jsc=1; j_add=0; jec=1
308 i1=1; i2=1; j1=1; j2=1;
310 allocate(local_buf_r8_kind(i1:i2,j1:j2,k_glob))
314 call mpp_scatter(isc+i_add, iec+i_add, jsc+j_add, jec+j_add, k_glob, bc_info%pelist, &
315 local_buf_r8_kind(i1:i2,j1:j2,:), global_buf_r8_kind, fileobj%is_root)
318 if (fileobj%is_root .and. .not. bc_info%data_on_file_root)
then
319 deallocate(local_buf_r8_kind)
324 vdata(i1:i2,j1:j2,:) = local_buf_r8_kind(i1:i2,j1:j2,:)
326 deallocate(local_buf_r8_kind)
327 if (fileobj%is_root)
deallocate(global_buf_r8_kind)
330 call error(
"scatter_data_bc_3d: unsupported type. Currently only r8_kind and r4_kinds are supported")
integer function mpp_pe()
Returns processor ID.
subroutine scatter_data_bc_2d(fileobj, varname, vdata, bc_info, unlim_dim_level, ignore_checksum)
subroutine scatter_data_bc_3d(fileobj, varname, vdata, bc_info, unlim_dim_level, ignore_checksum)