Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_geoamv / da_transform_xtoy_geoamv_adj.inc
blob52d0c4ae2d32bf0c913022e29161831a441feb31
1 subroutine da_transform_xtoy_geoamv_adj (grid, iv, jo_grad_y, jo_grad_x)
3    !-------------------------------------------------------------------------
4    ! Purpose: X to Y Transpose operator for Geo. AMVs
5    !    Updated for Analysis on Arakawa-C grid
6    !    Author: Syed RH Rizvi,  MMM/ESSL/NCAR,  Date: 10/22/2008
7    !-------------------------------------------------------------------------
9    implicit none
10    type(domain),  intent(in)    :: grid
11    type (iv_type), intent(in)    :: iv          ! obs. inc vector (o-b).
12    type (y_type) , intent(in)    :: jo_grad_y   ! grad_y(jo)
13    type (x_type) , intent(inout) :: jo_grad_x   ! grad_x(jo)
15    integer           :: n, k     ! Loop counter.
16    real, allocatable :: u(:,:)
17    real, allocatable :: v(:,:)
18    real, allocatable :: ub(:,:)
19    real, allocatable :: vb(:,:)
21    if (trace_use_dull) call da_trace_entry("da_transform_xtoy_geoamv_adj")
23    allocate (u(iv%info(geoamv)%max_lev,iv%info(geoamv)%n1:iv%info(geoamv)%n2))
24    allocate (v(iv%info(geoamv)%max_lev,iv%info(geoamv)%n1:iv%info(geoamv)%n2))
26    allocate (ub(iv%info(geoamv)%max_lev,iv%info(geoamv)%n1:iv%info(geoamv)%n2))
27    allocate (vb(iv%info(geoamv)%max_lev,iv%info(geoamv)%n1:iv%info(geoamv)%n2))
29    call da_interp_lin_3d (grid%xb%u, iv%info(geoamv), ub)
30    call da_interp_lin_3d (grid%xb%v, iv%info(geoamv), vb)
32    do n=iv%info(geoamv)%n1,iv%info(geoamv)%n2
33       do k = 1, iv%info(geoamv)%levels(n)
34           if(wind_sd_geoamv) then
35              call da_uv_to_sd_adj(jo_grad_y%geoamv(n)%u(k), &
36                                   jo_grad_y%geoamv(n)%v(k), u(k,n), v(k,n), ub(k,n), vb(k,n))
37           else
38              u(k,n) = jo_grad_y%geoamv(n)%u(k)
39              v(k,n) = jo_grad_y%geoamv(n)%v(k)
40           end if
41       end do
42    end do
44 #ifdef A2C
45    call da_interp_lin_3d_adj (jo_grad_x%u, iv%info(geoamv), u,'u')
46    call da_interp_lin_3d_adj (jo_grad_x%v, iv%info(geoamv), v,'v')
47 #else
48    call da_interp_lin_3d_adj (jo_grad_x%u, iv%info(geoamv), u)
49    call da_interp_lin_3d_adj (jo_grad_x%v, iv%info(geoamv), v)
50 #endif
52    deallocate (u)
53    deallocate (v)
54    deallocate (ub)
55    deallocate (vb)
57    if (trace_use_dull) call da_trace_exit("da_transform_xtoy_geoamv_adj")
59 end subroutine da_transform_xtoy_geoamv_adj