24 module fms_diag_input_buffer_mod
26 use platform_mod,
only: r8_kind, r4_kind, i4_kind, i8_kind
35 module procedure append_data_buffer_r4, append_data_buffer_r8
40 module procedure sum_data_buffer_r4, sum_data_buffer_r8
47 logical :: initialized
48 class(*),
allocatable :: buffer(:,:,:,:)
49 integer,
allocatable :: counter(:,:,:,:)
50 real(kind=r8_kind) :: weight
75 class(*),
pointer :: buffer(:,:,:,:)
86 real(kind=r8_kind),
pointer :: weight
96 class(*),
intent(in) :: input_data(:,:,:,:)
97 integer,
target,
intent(in) :: axis_ids(:)
99 character(len=128) :: err_msg
102 integer,
parameter :: ndims = 4
103 integer :: length(ndims)
105 integer,
pointer :: axis_id
113 naxes =
size(axis_ids)
114 axis_loop:
do a = 1,naxes
115 axis_id => axis_ids(a)
116 select type (axis => diag_axis(axis_id)%axis)
118 length(a) = axis%axis_length()
122 select type (input_data)
123 type is (real(r4_kind))
124 allocate(real(kind=r4_kind) :: this%buffer(length(1), length(2), length(3), length(4)))
125 this%buffer = 0.0_r4_kind
126 type is (real(r8_kind))
127 allocate(real(kind=r8_kind) :: this%buffer(length(1), length(2), length(3), length(4)))
128 this%buffer = 0.0_r8_kind
129 type is (
integer(i4_kind))
130 allocate(
integer(kind=i4_kind) :: this%buffer(length(1), length(2), length(3), length(4)))
131 this%buffer = 0_i4_kind
132 type is (
integer(i8_kind))
133 allocate(
integer(kind=i4_kind) :: this%buffer(length(1), length(2), length(3), length(4)))
134 this%buffer = 0_i8_kind
136 err_msg =
"The data input is not one of the supported types. &
137 &Only r4, r8, i4, and i8 types are supported."
140 this%weight = 1.0_r8_kind
141 this%initialized = .true.
142 allocate(this%counter(length(1), length(2), length(3), length(4)))
150 select type(buffer=>this%buffer)
151 type is (real(kind=r8_kind))
153 type is (real(kind=r4_kind))
164 this%send_data_time = time
174 rslt = this%send_data_time
179 function update_input_buffer_object(this, input_data, is, js, ks, ie, je, ke, mask_in, mask_out, &
180 mask_variant, var_is_masked) &
184 class(*),
intent(in) :: input_data(:,:,:,:)
185 integer,
intent(in) :: is, js, ks
186 integer,
intent(in) :: ie, je, ke
187 logical,
intent(in) :: mask_in(:,:,:,:)
188 logical,
intent(inout) :: mask_out(:,:,:,:)
189 logical,
intent(in) :: mask_variant
190 logical,
intent(in) :: var_is_masked
192 character(len=128) :: err_msg
194 if (mask_variant)
then
196 this%buffer(is:ie,js:je,ks:ke,:), input_data)
198 mask_out(is:ie,js:je,ks:ke,:) = mask_in
200 this%counter(is:ie,js:je,ks:ke,:), &
210 character(len=*),
intent(in) :: field_info
212 select type (input_data => this%buffer)
213 type is (real(kind=r4_kind))
214 input_data = input_data / this%counter(1,1,1,1)
215 type is (real(kind=r8_kind))
216 input_data = input_data / this%counter(1,1,1,1)
218 call mpp_error(fatal,
"prepare_input_buffer_object::"//trim(field_info)//&
219 " has only been implemented for real variables. Contact developers.")
228 logical,
intent(in) :: mask(:,:,:,:)
229 class(*),
intent(inout) :: data_out(:,:,:,:)
230 class(*),
intent(in) :: data_in(:,:,:,:)
231 integer,
intent(inout) :: counter(:,:,:,:)
232 logical,
intent(in) :: var_is_masked
234 character(len=128) :: err_msg
237 select type(data_out)
238 type is (real(kind=r8_kind))
239 select type (data_in)
240 type is (real(kind=r8_kind))
243 type is (real(kind=r4_kind))
244 select type (data_in)
245 type is (real(kind=r4_kind))
249 err_msg =
"sum_data_buffer_wrapper:: has only been implemented for real. Contact developers"
257 logical,
intent(inout) :: mask_out(:,:,:,:)
258 logical,
intent(in) :: mask_in(:,:,:,:)
259 class(*),
intent(inout) :: data_out(:,:,:,:)
260 class(*),
intent(in) :: data_in(:,:,:,:)
262 character(len=128) :: err_msg
265 select type(data_out)
266 type is (real(kind=r8_kind))
267 select type (data_in)
268 type is (real(kind=r8_kind))
271 type is (real(kind=r4_kind))
272 select type (data_in)
273 type is (real(kind=r4_kind))
277 err_msg =
"append_data_buffer:: has only been implemented for real. Contact developers"
287 class(*),
intent(in) :: input_data(:,:,:,:)
288 real(kind=r8_kind),
intent(in) :: weight
289 integer,
intent(in) :: is, js, ks
290 integer,
intent(in) :: ie, je, ke
292 character(len=128) :: err_msg
295 if (.not. this%initialized)
then
296 err_msg =
"The data buffer was never initiliazed. This shouldn't happen."
302 select type (input_data)
303 type is (real(kind=r4_kind))
304 select type (db => this%buffer)
305 type is (real(kind=r4_kind))
306 db(is:ie, js:je, ks:ke, :) = input_data
308 err_msg =
"The data buffer was not allocated to the correct type (r4_kind). This shouldn't happen"
311 type is (real(kind=r8_kind))
312 select type (db => this%buffer)
313 type is (real(kind=r8_kind))
314 db(is:ie, js:je, ks:ke, :) = input_data
316 err_msg =
"The data buffer was not allocated to the correct type (r8_kind). This shouldn't happen"
319 type is (
integer(kind=i4_kind))
320 select type (db => this%buffer)
321 type is (
integer(kind=i4_kind))
322 db(is:ie, js:je, ks:ke, :) = input_data
324 err_msg =
"The data buffer was not allocated to the correct type (i4_kind). This shouldn't happen"
327 type is (
integer(kind=i8_kind))
328 select type (db => this%buffer)
329 type is (
integer(kind=i8_kind))
330 db(is:ie, js:je, ks:ke, :) = input_data
332 err_msg =
"The data buffer was not allocated to the correct type (i8_kind). This shouldn't happen"
344 if (this%initialized)
then
351 #include "fms_diag_input_buffer_r4.fh"
352 #include "fms_diag_input_buffer_r8.fh"
355 end module fms_diag_input_buffer_mod
Type to represent amounts of time. Implemented as seconds and days to allow for larger intervals.
Type to hold the diag_axis (either subaxis or a full axis)
Type to hold the diagnostic axis description.