Update version info for release v4.6.1 (#2122)
[WRF.git] / var / da / da_airep / da_transform_xtoy_airep.inc
blob0a2d67b805ab408bb2bcdbe37012dc852a44ac94
1 subroutine da_transform_xtoy_airep (grid, iv, y)
3    !--------------------------------------------------------------------------
4    ! Purpose: TBD
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%grid%xa) (linear)
15    integer :: n,k        ! Loop counter.
17    real, allocatable :: u(:,:)
18    real, allocatable :: v(:,:)
19    real, allocatable :: t(:,:)
20    real, allocatable :: q(:,:)
22    real, allocatable :: ub(:,:)
23    real, allocatable :: vb(:,:)
25    if (trace_use_dull) call da_trace_entry("da_transform_xtoy_airep")
27    allocate (u(iv%info(airep)%max_lev,iv%info(airep)%n1:iv%info(airep)%n2))
28    allocate (v(iv%info(airep)%max_lev,iv%info(airep)%n1:iv%info(airep)%n2))
29    allocate (t(iv%info(airep)%max_lev,iv%info(airep)%n1:iv%info(airep)%n2))
30    allocate (q(iv%info(airep)%max_lev,iv%info(airep)%n1:iv%info(airep)%n2))
31   
32    allocate (ub(iv%info(airep)%max_lev,iv%info(airep)%n1:iv%info(airep)%n2))
33    allocate (vb(iv%info(airep)%max_lev,iv%info(airep)%n1:iv%info(airep)%n2))
35 #ifdef A2C
36    call da_interp_lin_3d (grid%xa%u, iv%info(airep), u,'u')
37    call da_interp_lin_3d (grid%xa%v, iv%info(airep), v,'v')
38 #else
39    call da_interp_lin_3d (grid%xa%u, iv%info(airep), u)
40    call da_interp_lin_3d (grid%xa%v, iv%info(airep), v)
41 #endif
42    call da_interp_lin_3d (grid%xa%t, iv%info(airep), t)
43    call da_interp_lin_3d (grid%xa%q, iv%info(airep), q)
45    call da_interp_lin_3d (grid%xb%u, iv%info(airep), ub)
46    call da_interp_lin_3d (grid%xb%v, iv%info(airep), vb)
48    do n=iv%info(airep)%n1,iv%info(airep)%n2
49       do k = 1, iv%info(airep)%levels(n)
50          if(wind_sd_airep) then
51              call da_uv_to_sd_lin(y%airep(n)%u(k),y%airep(n)%v(k),u(k,n),v(k,n),ub(k,n),vb(k,n))
52          else
53              y%airep(n)%u(k) = u(k,n)
54              y%airep(n)%v(k) = v(k,n)
55          endif
56          y%airep(n)%t(:) = t(1:size(y%airep(n)%t),n)
57          y%airep(n)%q(:) = q(1:size(y%airep(n)%q),n)
58       end do
59    end do
60   
61    deallocate (u)
62    deallocate (v)
63    deallocate (t)
64    deallocate (q)
65    deallocate (ub)
66    deallocate (vb)
68    if (trace_use_dull) call da_trace_exit("da_transform_xtoy_airep")
70 end subroutine da_transform_xtoy_airep