updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_airsr / da_transform_xtoy_airsr_adj.inc
blob796c4bb7d7f562450feae00a46325b3ad2299359
1 subroutine da_transform_xtoy_airsr_adj(iv, jo_grad_y, jo_grad_x)
3    !----------------------------------------------------------------------
4    ! Purpose: Does adjoint computation at AIRS retrieval locations
5    !----------------------------------------------------------------------
7    implicit none
9    type (iv_type), intent(in)    :: iv          ! obs. inc vector (o-b).
10    type (y_type) , intent(in)    :: jo_grad_y   ! grad_y(jo)
11    type (x_type) , intent(inout) :: jo_grad_x   ! grad_x(jo)
13    integer :: n  ! Loop counter.
15    real, allocatable :: t(:,:)
16    real, allocatable :: q(:,:)           
18    if (trace_use_dull) call da_trace_entry("da_transform_xtoy_airsr_adj")
20    allocate (t(iv%info(airsr)%max_lev,iv%info(airsr)%n1:iv%info(airsr)%n2))
21    allocate (q(iv%info(airsr)%max_lev,iv%info(airsr)%n1:iv%info(airsr)%n2))
23    do n=iv%info(airsr)%n1,iv%info(airsr)%n2
24       t(1:size(jo_grad_y%airsr(n)%t),n)  = jo_grad_y%airsr(n)%t
25       q(1:size(jo_grad_y%airsr(n)%q),n)  = jo_grad_y%airsr(n)%q
26    end do
28    ! [1.1] Adjoint feedback from Y to X for u and v:
30    call da_interp_lin_3d_adj (jo_grad_x%t, iv%info(airsr), t)
31    call da_interp_lin_3d_adj (jo_grad_x%q, iv%info(airsr), q)
33    deallocate (t)
34    deallocate (q)
36    if (trace_use_dull) call da_trace_exit("da_transform_xtoy_airsr_adj")
38 end subroutine da_transform_xtoy_airsr_adj