da_transform_xtoy_sound_adj.inc

References to this file elsewhere.
1 subroutine da_transform_xtoy_sound_adj(iv, jo_grad_y, jo_grad_x)
2 
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
6 
7    implicit none
8 
9    type (iv_type), intent(in)    :: iv          ! obs. inc vector (o-b).
10    type (y_type) , intent(in)    :: jo_grad_y   ! grad_y(jo)
11    type (x_type) , intent(inout) :: jo_grad_x   ! grad_x(jo)
12 
13    integer :: n
14 
15    real, allocatable :: u(:,:)
16    real, allocatable :: v(:,:)
17    real, allocatable :: t(:,:)
18    real, allocatable :: q(:,:)
19 
20    if (trace_use_dull) call da_trace_entry("da_transform_xtoy_sound_adj")
21 
22    allocate (u(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%n2))
23    allocate (v(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%n2))
24    allocate (t(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%n2))
25    allocate (q(iv%info(sound)%max_lev,iv%info(sound)%n1:iv%info(sound)%n2))
26 
27    do n=iv%info(sound)%n1,iv%info(sound)%n2
28       u(1:size(jo_grad_y%sound(n)%u),n) = jo_grad_y%sound(n)%u(:)
29       v(1:size(jo_grad_y%sound(n)%v),n) = jo_grad_y%sound(n)%v(:)
30       t(1:size(jo_grad_y%sound(n)%t),n) = jo_grad_y%sound(n)%t(:)
31       q(1:size(jo_grad_y%sound(n)%q),n) = jo_grad_y%sound(n)%q(:)
32    end do
33 
34    call da_interp_lin_3d_adj (jo_grad_x%u, iv%info(sound), u)
35    call da_interp_lin_3d_adj (jo_grad_x%v, iv%info(sound), v)
36    call da_interp_lin_3d_adj (jo_grad_x%t, iv%info(sound), t)
37    call da_interp_lin_3d_adj (jo_grad_x%q, iv%info(sound), q)
38 
39    deallocate (u)
40    deallocate (v)
41    deallocate (t)
42    deallocate (q)
43 
44    if (trace_use_dull) call da_trace_exit("da_transform_xtoy_sound_adj")
45 
46 end subroutine da_transform_xtoy_sound_adj
47 
48