da_transform_xtoy_airep_adj.inc

References to this file elsewhere.
1 subroutine da_transform_xtoy_airep_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 
19    if (trace_use_dull) call da_trace_entry("da_transform_xtoy_airep_adj")
20 
21    allocate (u(iv%info(airep)%max_lev,iv%info(airep)%n1:iv%info(airep)%n2))
22    allocate (v(iv%info(airep)%max_lev,iv%info(airep)%n1:iv%info(airep)%n2))
23    allocate (t(iv%info(airep)%max_lev,iv%info(airep)%n1:iv%info(airep)%n2))
24 
25    do n=iv%info(airep)%n1,iv%info(airep)%n2
26       u(1:size(jo_grad_y%airep(n)%u),n) = jo_grad_y%airep(n)%u(:)
27       v(1:size(jo_grad_y%airep(n)%v),n) = jo_grad_y%airep(n)%v(:)
28       t(1:size(jo_grad_y%airep(n)%t),n) = jo_grad_y%airep(n)%t(:)
29    end do
30 
31    call da_interp_lin_3d_adj(jo_grad_x%u, iv%info(airep), u)
32    call da_interp_lin_3d_adj(jo_grad_x%v, iv%info(airep), v)
33    call da_interp_lin_3d_adj(jo_grad_x%t, iv%info(airep), t)
34 
35    deallocate (u)
36    deallocate (v)
37    deallocate (t)
38 
39    if (trace_use_dull) call da_trace_exit("da_transform_xtoy_airep_adj")
40 
41 end subroutine da_transform_xtoy_airep_adj
42 
43