Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_tamdar / da_transform_xtoy_tamdar.inc
blob2dd286a79d89451c536076e602cfd32fa302c6af
1 subroutine da_transform_xtoy_tamdar (grid, iv, y)
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
7    implicit none
9    type (domain),     intent(in)    :: grid
10    type (iv_type),    intent(in)    :: iv       ! Innovation vector (O-B).
11    type (y_type),     intent(inout) :: y        ! y = h (grid%xa) (linear)
13    real, allocatable :: u(:,:)
14    real, allocatable :: v(:,:)
15    real, allocatable :: t(:,:)
16    real, allocatable :: q(:,:)
18    real, allocatable :: ub(:,:)
19    real, allocatable :: vb(:,:)
21    integer :: n, k
23    if (trace_use_dull) call da_trace_entry("da_transform_xtoy_tamdar")
25    allocate (u(iv%info(tamdar)%max_lev,iv%info(tamdar)%n1:iv%info(tamdar)%n2))
26    allocate (v(iv%info(tamdar)%max_lev,iv%info(tamdar)%n1:iv%info(tamdar)%n2))
27    allocate (t(iv%info(tamdar)%max_lev,iv%info(tamdar)%n1:iv%info(tamdar)%n2))
28    allocate (q(iv%info(tamdar)%max_lev,iv%info(tamdar)%n1:iv%info(tamdar)%n2))
30    allocate (ub(iv%info(tamdar)%max_lev,iv%info(tamdar)%n1:iv%info(tamdar)%n2))
31    allocate (vb(iv%info(tamdar)%max_lev,iv%info(tamdar)%n1:iv%info(tamdar)%n2))
33 #ifdef A2C  
34    call da_interp_lin_3d (grid%xa%u, iv%info(tamdar), u, 'u')
35    call da_interp_lin_3d (grid%xa%v, iv%info(tamdar), v, 'v')
36 #else
37    call da_interp_lin_3d (grid%xa%u, iv%info(tamdar), u)
38    call da_interp_lin_3d (grid%xa%v, iv%info(tamdar), v)
39 #endif
40    call da_interp_lin_3d (grid%xa%t, iv%info(tamdar), t)
41    call da_interp_lin_3d (grid%xa%q, iv%info(tamdar), q)
43    call da_interp_lin_3d (grid%xb%u, iv%info(tamdar), ub)
44    call da_interp_lin_3d (grid%xb%v, iv%info(tamdar), vb)
46    do n=iv%info(tamdar)%n1,iv%info(tamdar)%n2
47       do k = 1, iv%info(tamdar)%levels(n)
48          if(wind_sd_tamdar) then
49             call da_uv_to_sd_lin(y%tamdar(n)%u(k),y%tamdar(n)%v(k),u(k,n),v(k,n),ub(k,n),vb(k,n))
50          else
51             y%tamdar(n)%u(k) = u(k,n)
52             y%tamdar(n)%v(k) = v(k,n)
53          end if
54       end do
55       y%tamdar(n)%t(:) = t(1:size(y%tamdar(n)%t),n)
56       y%tamdar(n)%q(:) = q(1:size(y%tamdar(n)%q),n)
57    end do
59    deallocate (u)
60    deallocate (v)
61    deallocate (t)
62    deallocate (q)
63    deallocate (ub)
64    deallocate (vb)
66    if (trace_use_dull) call da_trace_exit("da_transform_xtoy_tamdar")
68 end subroutine da_transform_xtoy_tamdar