Skip to content

Commit

Permalink
Allow setting of the fill value when expanding the buffer
Browse files Browse the repository at this point in the history
  • Loading branch information
ashao committed Jan 23, 2025
1 parent cc27731 commit c1ab465
Show file tree
Hide file tree
Showing 2 changed files with 14 additions and 3 deletions.
13 changes: 11 additions & 2 deletions src/framework/MOM_diag_buffers.F90
Original file line number Diff line number Diff line change
Expand Up @@ -16,13 +16,15 @@ module MOM_diag_buffers
integer :: js !< The start slot of the array j-direction
integer :: ie !< The end slot of the array i-direction
integer :: je !< The end slot of the array j-direction
real :: fill_value = 0. !< Set the fill value to use when growing the buffer

integer, allocatable, dimension(:) :: ids !< List of diagnostic ids whose slot corresponds to the row in the buffer
integer :: length = 0 !< The number of slots in the buffer

contains

procedure(a_grow), deferred :: grow !< Increase the size of the buffer
procedure, public :: set_fill_value !< Set the fill value to use when growing the buffer
procedure, public :: check_capacity_by_id !< Check the size size of the buffer and increase if necessary
procedure, public :: set_horizontal_extents !< Define the horizontal extents of the arrays
procedure, public :: mark_available !< Mark that a slot in the buffer can be reused
Expand Down Expand Up @@ -60,6 +62,13 @@ subroutine a_grow(this)
class(diag_buffer_base), intent(inout) :: this !< The diagnostic buffer
end subroutine

!> Set the fill value to use when growing the buffer
subroutine set_fill_value(this, fill_value)
class(diag_buffer_base), intent(inout) :: this !< The diagnostic buffer
real, intent(in) :: fill_value !< The fill value to use when growing the buffer [arbitrary]

end subroutine set_fill_value

!> Mark a slot in the buffer as unused based on a diagnostic id. For example,
!! the data in that slot has already been consumed and can thus be overwritten
subroutine mark_available(this, id)
Expand Down Expand Up @@ -154,7 +163,7 @@ subroutine grow_2d(this)
call this%grow_ids()

n = this%length
allocate(temp(n+1, this%is:this%ie, this%js:this%je), source=0.)
allocate(temp(n+1, this%is:this%ie, this%js:this%je), source=this%fill_value)
if (n>0) temp(1:n,:,:) = this%buffer(:,:,:)
call move_alloc(temp, this%buffer)
this%length = this%length + 1
Expand Down Expand Up @@ -183,7 +192,7 @@ subroutine grow_3d(this)
call this%grow_ids()

n = this%length
allocate(temp(n+1, this%is:this%ie, this%js:this%je, this%ks:this%ke), source=0.)
allocate(temp(n+1, this%is:this%ie, this%js:this%je, this%ks:this%ke), source=this%fill_value)
if (n>0) temp(1:n,:,:,:) = this%buffer(:,:,:,:)
call move_alloc(temp, this%buffer)
this%length = this%length + 1
Expand Down
4 changes: 3 additions & 1 deletion src/framework/MOM_diag_mediator.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1109,6 +1109,7 @@ subroutine define_axes_group(diag_cs, handles, axes, nz, vertical_coordinate_num
if (axes%is_q_point) axes%mask2d => diag_cs%mask2dBu
call axes%buffer_2d%set_horizontal_extents(lbound(axes%mask2d,1), ubound(axes%mask2d,1), &
lbound(axes%mask2d,2), ubound(axes%mask2d,2))
call axes%buffer_2d%set_fill_value(diag_cs%missing_value)
endif
! A static 3d mask for non-native coordinates can only be setup when a grid is available
axes%mask3d => null()
Expand All @@ -1128,6 +1129,7 @@ subroutine define_axes_group(diag_cs, handles, axes, nz, vertical_coordinate_num
call axes%buffer_3d%set_horizontal_extents(is=lbound(axes%mask3d,1), ie=ubound(axes%mask3d,1), &
js=lbound(axes%mask3d,2), je=ubound(axes%mask3d,2))
call axes%buffer_3d%set_vertical_extent(ks=lbound(axes%mask3d,3), ke=ubound(axes%mask3d,3))
call axes%buffer_3d%set_fill_value(diag_cs%missing_value)
endif


Expand Down Expand Up @@ -1931,7 +1933,7 @@ subroutine post_data_3d_by_point(diag_field_id, field, diag_cs, i, j, k)

diag => diag_cs%diags(diag_field_id)
buffer_slot = diag%axes%buffer_3d%find_buffer_slot(diag_field_id)
diag%axes%buffer_3d%buffer(buffer_slot, i, j, k) = field
diag%axes%buffer_3d%buffer(buffer_slot,i,j,k) = field
end subroutine post_data_3d_by_point

!> Post the final buffer using the standard post_data interface
Expand Down

0 comments on commit c1ab465

Please sign in to comment.