updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_sound / da_transform_xtoy_sound_adj.inc
blobb2ef663dfbc8c115bdf51750a72b4ac815ab6730
1 subroutine da_transform_xtoy_sound_adj(grid, iv, jo_grad_y, jo_grad_x)
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
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
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_sound_adj")
27    allocate (u(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%n2))
28    allocate (v(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%n2))
29    allocate (t(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%n2))
30    allocate (q(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%n2))
32    allocate (ub(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%n2))
33    allocate (vb(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%n2))
35    call da_interp_lin_3d (grid%xb%u, iv%info(sound), ub)
36    call da_interp_lin_3d (grid%xb%v, iv%info(sound), vb)
38    do n=iv%info(sound)%n1,iv%info(sound)%n2
39        do k = 1, iv%info(sound)%levels(n)
40          if(wind_sd_sound) then
41             call da_uv_to_sd_adj(jo_grad_y%sound(n)%u(k), &
42                                  jo_grad_y%sound(n)%v(k), u(k,n), v(k,n), ub(k,n), vb(k,n))
43          else
44             u(k,n) = jo_grad_y%sound(n)%u(k)
45             v(k,n) = jo_grad_y%sound(n)%v(k)
46          end if
47       end do
48       t(1:size(jo_grad_y%sound(n)%t),n) = jo_grad_y%sound(n)%t(:)
49       q(1:size(jo_grad_y%sound(n)%q),n) = jo_grad_y%sound(n)%q(:)
50    end do
52 #ifdef A2C
53    call da_interp_lin_3d_adj (jo_grad_x%u, iv%info(sound), u,'u')
54    call da_interp_lin_3d_adj (jo_grad_x%v, iv%info(sound), v,'v')
55 #else
56    call da_interp_lin_3d_adj (jo_grad_x%u, iv%info(sound), u)
57    call da_interp_lin_3d_adj (jo_grad_x%v, iv%info(sound), v)
58 #endif
59    call da_interp_lin_3d_adj (jo_grad_x%t, iv%info(sound), t)
60    call da_interp_lin_3d_adj (jo_grad_x%q, iv%info(sound), q)
61    
62    deallocate (u)
63    deallocate (v)
64    deallocate (t)
65    deallocate (q)
66    deallocate (ub)
67    deallocate (vb)
69    if (trace_use_dull) call da_trace_exit("da_transform_xtoy_sound_adj")
71 end subroutine da_transform_xtoy_sound_adj