Update version info for release v4.6.1 (#2122)
[WRF.git] / var / da / da_pseudo / da_transform_xtoy_pseudo_adj.inc
blob2e0990c8c245e398e50184e85fcc420e7e173bfe
1 subroutine da_transform_xtoy_pseudo_adj(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
11    type (iv_type), intent(in)    :: iv          ! obs. inc vector (o-b).
12    type (y_type) , intent(inout) :: jo_grad_y   ! grad_y(jo)
13    type (x_type) , intent(inout) :: jo_grad_x   ! grad_x(jo)
15    integer :: n        ! Loop counter.
17    real, allocatable :: u(:,:)
18    real, allocatable :: v(:,:)
19    real, allocatable :: q(:,:)
20    real, allocatable :: p(:,:)
21    real, allocatable :: t(:,:)
22    real, allocatable :: qw(:,:)
24    if (trace_use_dull) call da_trace_entry("da_transform_xtoy_pseudo_adj")
26    allocate (u(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
27    allocate (v(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
28    allocate (q(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
29    allocate (p(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
30    allocate (t(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
31    allocate (qw(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
33    do n=iv%info(pseudo)%n1,iv%info(pseudo)%n2
34       u(1,n) = jo_grad_y%pseudo(n)%u
35       v(1,n) = jo_grad_y%pseudo(n)%v
37       if ( len_trim(adjustl(pseudo_var)) == 1 ) then
38          q(1,n) = jo_grad_y%pseudo(n)%q
39       else
40          select case(pseudo_var(1:3))
41          case ('qcw', 'QCW','qci', 'QCI','qrn', 'QRN', 'qsn', 'QSN', 'qgr', 'QGR')
42             qw(1,n) = jo_grad_y%pseudo(n)%q
43          end select
44       end if
46       p(1,n) = jo_grad_y%pseudo(n)%p
47       t(1,n) = jo_grad_y%pseudo(n)%t
48    end do
50 #ifdef A2C
51    call da_interp_lin_3d_adj(jo_grad_x%u, iv%info(pseudo), u,'u')
52    call da_interp_lin_3d_adj(jo_grad_x%v, iv%info(pseudo), v,'v')
53 #else
54    call da_interp_lin_3d_adj(jo_grad_x%u, iv%info(pseudo), u)
55    call da_interp_lin_3d_adj(jo_grad_x%v, iv%info(pseudo), v)
56 #endif
57    call da_interp_lin_3d_adj(jo_grad_x%q, iv%info(pseudo), q)
58    call da_interp_lin_3d_adj(jo_grad_x%p, iv%info(pseudo), p)
59    call da_interp_lin_3d_adj(jo_grad_x%t, iv%info(pseudo), t)
61    select case(pseudo_var(1:3))
62    case ('qcw', 'QCW')
63       call da_interp_lin_3d_adj(jo_grad_x%qcw, iv%info(pseudo), qw)
64    case ('qci', 'QCI')
65       call da_interp_lin_3d_adj(jo_grad_x%qci, iv%info(pseudo), qw)
66    case ('qrn', 'QRN')
67       call da_interp_lin_3d_adj(jo_grad_x%qrn, iv%info(pseudo), qw)
68    case ('qsn', 'QSN')
69       call da_interp_lin_3d_adj(jo_grad_x%qsn, iv%info(pseudo), qw)
70    case ('qgr', 'QGR')
71       call da_interp_lin_3d_adj(jo_grad_x%qgr, iv%info(pseudo), qw)
72    end select
74    deallocate (u)
75    deallocate (v)
76    deallocate (q)
77    deallocate (p)
78    deallocate (t)
79    deallocate (qw)
81    if (trace_use_dull) call da_trace_exit("da_transform_xtoy_pseudo_adj")
83 end subroutine da_transform_xtoy_pseudo_adj