Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_geoamv / da_transform_xtoy_geoamv.inc
blobfa1b339226f013a4f5487e6827365aa0da9ed86c
1 subroutine da_transform_xtoy_geoamv (grid, iv, y)
3    !-------------------------------------------------------------------------
4    ! Purpose: X to Y Transform 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
11    type (domain),  intent(in)    :: grid
12    type (iv_type), intent(in)    :: iv       ! Innovation vector (O-B).
13    type (y_type),  intent(inout) :: y        ! y = h (grid%xa)
15    integer           :: n,k
16    real, allocatable :: u(:,:)
17    real, allocatable :: v(:,:)
19    real, allocatable :: ub(:,:)
20    real, allocatable :: vb(:,:)
22    if (trace_use_dull) call da_trace_entry("da_transform_xtoy_geoamv")
24    allocate (u(iv%info(geoamv)%max_lev,iv%info(geoamv)%n1:iv%info(geoamv)%n2))
25    allocate (v(iv%info(geoamv)%max_lev,iv%info(geoamv)%n1:iv%info(geoamv)%n2))
27    allocate (ub(iv%info(geoamv)%max_lev,iv%info(geoamv)%n1:iv%info(geoamv)%n2))
28    allocate (vb(iv%info(geoamv)%max_lev,iv%info(geoamv)%n1:iv%info(geoamv)%n2))
30 #ifdef A2C
31    call da_interp_lin_3d (grid%xa%u, iv%info(geoamv), u,'u')
32    call da_interp_lin_3d (grid%xa%v, iv%info(geoamv), v,'v')
33 #else
34    call da_interp_lin_3d (grid%xa%u, iv%info(geoamv), u)
35    call da_interp_lin_3d (grid%xa%v, iv%info(geoamv), v)
36 #endif
38    call da_interp_lin_3d (grid%xb%u, iv%info(geoamv), ub)
39    call da_interp_lin_3d (grid%xb%v, iv%info(geoamv), vb)
41    do n=iv%info(geoamv)%n1,iv%info(geoamv)%n2
42       do k = 1, iv%info(geoamv)%levels(n)
43          if(wind_sd_geoamv) then
44             call da_uv_to_sd_lin(y%geoamv(n)%u(k),y%geoamv(n)%v(k),u(k,n),v(k,n),ub(k,n),vb(k,n))
45          else
46             y%geoamv(n)%u(k) = u(k,n)
47             y%geoamv(n)%v(k) = v(k,n)
48          end if
49        end do
50    end do
52    deallocate (u)
53    deallocate (v)
54    deallocate (vb)
55    deallocate (ub)
57    if (trace_use_dull) call da_trace_exit("da_transform_xtoy_geoamv")
59 end subroutine da_transform_xtoy_geoamv