23 module fms_diag_input_buffer_mod
25 use platform_mod,
only: r8_kind, r4_kind, i4_kind, i8_kind
34 module procedure append_data_buffer_r4, append_data_buffer_r8
39 module procedure sum_data_buffer_r4, sum_data_buffer_r8
46 logical :: initialized
47 class(*),
allocatable :: buffer(:,:,:,:)
48 integer,
allocatable :: counter(:,:,:,:)
49 real(kind=r8_kind) :: weight
74 class(*),
pointer :: buffer(:,:,:,:)
85 real(kind=r8_kind),
pointer :: weight
95 class(*),
intent(in) :: input_data(:,:,:,:)
96 integer,
target,
intent(in) :: axis_ids(:)
98 character(len=128) :: err_msg
101 integer,
parameter :: ndims = 4
102 integer :: length(ndims)
104 integer,
pointer :: axis_id
112 naxes =
size(axis_ids)
113 axis_loop:
do a = 1,naxes
114 axis_id => axis_ids(a)
115 select type (axis => diag_axis(axis_id)%axis)
117 length(a) = axis%axis_length()
121 select type (input_data)
122 type is (real(r4_kind))
123 allocate(real(kind=r4_kind) :: this%buffer(length(1), length(2), length(3), length(4)))
124 this%buffer = 0.0_r4_kind
125 type is (real(r8_kind))
126 allocate(real(kind=r8_kind) :: this%buffer(length(1), length(2), length(3), length(4)))
127 this%buffer = 0.0_r8_kind
128 type is (
integer(i4_kind))
129 allocate(
integer(kind=i4_kind) :: this%buffer(length(1), length(2), length(3), length(4)))
130 this%buffer = 0_i4_kind
131 type is (
integer(i8_kind))
132 allocate(
integer(kind=i4_kind) :: this%buffer(length(1), length(2), length(3), length(4)))
133 this%buffer = 0_i8_kind
135 err_msg =
"The data input is not one of the supported types. &
136 &Only r4, r8, i4, and i8 types are supported."
139 this%weight = 1.0_r8_kind
140 this%initialized = .true.
141 allocate(this%counter(length(1), length(2), length(3), length(4)))
149 select type(buffer=>this%buffer)
150 type is (real(kind=r8_kind))
152 type is (real(kind=r4_kind))
163 this%send_data_time = time
173 rslt = this%send_data_time
178 function update_input_buffer_object(this, input_data, is, js, ks, ie, je, ke, mask_in, mask_out, &
179 mask_variant, var_is_masked) &
183 class(*),
intent(in) :: input_data(:,:,:,:)
184 integer,
intent(in) :: is, js, ks
185 integer,
intent(in) :: ie, je, ke
186 logical,
intent(in) :: mask_in(:,:,:,:)
187 logical,
intent(inout) :: mask_out(:,:,:,:)
188 logical,
intent(in) :: mask_variant
189 logical,
intent(in) :: var_is_masked
191 character(len=128) :: err_msg
193 if (mask_variant)
then
195 this%buffer(is:ie,js:je,ks:ke,:), input_data)
197 mask_out(is:ie,js:je,ks:ke,:) = mask_in
199 this%counter(is:ie,js:je,ks:ke,:), &
209 character(len=*),
intent(in) :: field_info
211 select type (input_data => this%buffer)
212 type is (real(kind=r4_kind))
213 input_data = input_data / this%counter(1,1,1,1)
214 type is (real(kind=r8_kind))
215 input_data = input_data / this%counter(1,1,1,1)
217 call mpp_error(fatal,
"prepare_input_buffer_object::"//trim(field_info)//&
218 " has only been implemented for real variables. Contact developers.")
227 logical,
intent(in) :: mask(:,:,:,:)
228 class(*),
intent(inout) :: data_out(:,:,:,:)
229 class(*),
intent(in) :: data_in(:,:,:,:)
230 integer,
intent(inout) :: counter(:,:,:,:)
231 logical,
intent(in) :: var_is_masked
233 character(len=128) :: err_msg
236 select type(data_out)
237 type is (real(kind=r8_kind))
238 select type (data_in)
239 type is (real(kind=r8_kind))
242 type is (real(kind=r4_kind))
243 select type (data_in)
244 type is (real(kind=r4_kind))
248 err_msg =
"sum_data_buffer_wrapper:: has only been implemented for real. Contact developers"
256 logical,
intent(inout) :: mask_out(:,:,:,:)
257 logical,
intent(in) :: mask_in(:,:,:,:)
258 class(*),
intent(inout) :: data_out(:,:,:,:)
259 class(*),
intent(in) :: data_in(:,:,:,:)
261 character(len=128) :: err_msg
264 select type(data_out)
265 type is (real(kind=r8_kind))
266 select type (data_in)
267 type is (real(kind=r8_kind))
270 type is (real(kind=r4_kind))
271 select type (data_in)
272 type is (real(kind=r4_kind))
276 err_msg =
"append_data_buffer:: has only been implemented for real. Contact developers"
286 class(*),
intent(in) :: input_data(:,:,:,:)
287 real(kind=r8_kind),
intent(in) :: weight
288 integer,
intent(in) :: is, js, ks
289 integer,
intent(in) :: ie, je, ke
291 character(len=128) :: err_msg
294 if (.not. this%initialized)
then
295 err_msg =
"The data buffer was never initiliazed. This shouldn't happen."
301 select type (input_data)
302 type is (real(kind=r4_kind))
303 select type (db => this%buffer)
304 type is (real(kind=r4_kind))
305 db(is:ie, js:je, ks:ke, :) = input_data
307 err_msg =
"The data buffer was not allocated to the correct type (r4_kind). This shouldn't happen"
310 type is (real(kind=r8_kind))
311 select type (db => this%buffer)
312 type is (real(kind=r8_kind))
313 db(is:ie, js:je, ks:ke, :) = input_data
315 err_msg =
"The data buffer was not allocated to the correct type (r8_kind). This shouldn't happen"
318 type is (
integer(kind=i4_kind))
319 select type (db => this%buffer)
320 type is (
integer(kind=i4_kind))
321 db(is:ie, js:je, ks:ke, :) = input_data
323 err_msg =
"The data buffer was not allocated to the correct type (i4_kind). This shouldn't happen"
326 type is (
integer(kind=i8_kind))
327 select type (db => this%buffer)
328 type is (
integer(kind=i8_kind))
329 db(is:ie, js:je, ks:ke, :) = input_data
331 err_msg =
"The data buffer was not allocated to the correct type (i8_kind). This shouldn't happen"
343 if (this%initialized)
then
350 #include "fms_diag_input_buffer_r4.fh"
351 #include "fms_diag_input_buffer_r8.fh"
354 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.