Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

capture QC changes if posterior FOs skipped #430

Closed
wants to merge 1 commit into from
Closed
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
56 changes: 55 additions & 1 deletion assimilation_code/modules/assimilation/filter_mod.f90
Original file line number Diff line number Diff line change
Expand Up @@ -1051,6 +1051,10 @@ subroutine filter_main()

call trace_message('After posterior obs space diagnostics')
else
! call this alternate routine to collect any updated QC values that may
! have been set in the assimilation loop and copy them to the outgoing obs seq
call obs_space_sync_QCs(obs_fwd_op_ens_handle, seq, keys, num_obs_in_set, &
OBS_GLOBAL_QC_COPY, DART_qc_index)
call deallocate_single_copy(obs_fwd_op_ens_handle, prior_qc_copy)
endif

Expand Down Expand Up @@ -1663,7 +1667,7 @@ subroutine obs_space_diagnostics(obs_fwd_op_ens_handle, qc_ens_handle, ens_size,
OBS_MEAN_START, OBS_VAR_START, OBS_GLOBAL_QC_COPY, OBS_VAL_COPY, &
OBS_ERR_VAR_COPY, DART_qc_index, do_post)

! Do prior observation space diagnostics on the set of obs corresponding to keys
! Do observation space diagnostics on the set of obs corresponding to keys

type(ensemble_type), intent(inout) :: obs_fwd_op_ens_handle, qc_ens_handle
integer, intent(in) :: ens_size
Expand Down Expand Up @@ -1761,6 +1765,56 @@ end subroutine obs_space_diagnostics

!-------------------------------------------------------------------------

subroutine obs_space_sync_QCs(obs_fwd_op_ens_handle, &
seq, keys, num_obs_in_set, OBS_GLOBAL_QC_COPY, DART_qc_index)

! If QCs were updated in the assimilation loop but posterior forward operators
! are not being computed, collect any updated QCs into the output obs_seq.

type(ensemble_type), intent(inout) :: obs_fwd_op_ens_handle
integer, intent(in) :: num_obs_in_set
integer, intent(in) :: keys(num_obs_in_set)
type(obs_sequence_type), intent(inout) :: seq
integer, intent(in) :: OBS_GLOBAL_QC_COPY
integer, intent(in) :: DART_qc_index

integer :: j
integer :: io_task, my_task
real(r8), allocatable :: obs_temp(:)
real(r8) :: rvalue(1)

! this is a query routine to return which task has
! logical processing element 0 in this ensemble.
io_task = map_pe_to_task(obs_fwd_op_ens_handle, 0)
my_task = my_task_id()

! Make var complete for get_copy() calls below.
! Optimize: Could we use a gather instead of a transpose and get copy?
call all_copies_to_all_vars(obs_fwd_op_ens_handle)

! allocate temp space for sending data only on the task that will
! write the obs_seq.final file
if (my_task == io_task) then
allocate(obs_temp(num_obs_in_set))
else ! TJH: this change became necessary when using Intel 19.0.5 ...
allocate(obs_temp(1))
endif

! Update the qc global value
call get_copy(io_task, obs_fwd_op_ens_handle, OBS_GLOBAL_QC_COPY, obs_temp)
if(my_task == io_task) then
do j = 1, obs_fwd_op_ens_handle%num_vars
rvalue(1) = obs_temp(j)
call replace_qc(seq, keys(j), rvalue, DART_qc_index)
end do
endif

deallocate(obs_temp)

end subroutine obs_space_sync_QCs

!-------------------------------------------------------------------------

subroutine filter_sync_keys_time(ens_handle, key_bounds, num_obs_in_set, time1, time2)

integer, intent(inout) :: key_bounds(2), num_obs_in_set
Expand Down