Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_pseudo / da_transform_xtoy_pseudo.inc
blob7538f5ee5b3e40fa3d23b57606e72759d5feb0b0
1 subroutine da_transform_xtoy_pseudo(grid, iv, y)
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 (domain),  intent(inout) :: grid
12    type (iv_type), intent(in)    :: iv       ! Innovation vector (O-B).
13    type (y_type),  intent(inout) :: y        ! y = h (grid%xa)
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 !  working array for 5 hydrometeor variables
23    real, allocatable :: qw(:,:)
25    if (trace_use_dull) call da_trace_entry("da_transform_xtoy_pseudo")
27    allocate (u(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
28    allocate (v(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
29    allocate (q(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
30    allocate (p(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
31    allocate (t(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
32    allocate (qw(1,iv%info(pseudo)%n1:iv%info(pseudo)%n2))
34 #ifdef A2C
35    call da_interp_lin_3d(grid%xa%u, iv%info(pseudo), u,'u')
36    call da_interp_lin_3d(grid%xa%v, iv%info(pseudo), v,'v')
37 #else
38    call da_interp_lin_3d(grid%xa%u, iv%info(pseudo), u)
39    call da_interp_lin_3d(grid%xa%v, iv%info(pseudo), v)
40 #endif
41    call da_interp_lin_3d(grid%xa%q, iv%info(pseudo), q)
42    call da_interp_lin_3d(grid%xa%p, iv%info(pseudo), p)
43    call da_interp_lin_3d(grid%xa%t, iv%info(pseudo), t)
45    select case(pseudo_var(1:3))
46    case ('qcw', 'QCW')
47         call da_interp_lin_3d(grid%xa%qcw, iv%info(pseudo), qw)
48    case ('qci', 'QCI')
49         call da_interp_lin_3d(grid%xa%qci, iv%info(pseudo), qw)
50    case ('qrn', 'QRN')
51         call da_interp_lin_3d(grid%xa%qrn, iv%info(pseudo), qw)
52    case ('qsn', 'QSN')
53         call da_interp_lin_3d(grid%xa%qsn, iv%info(pseudo), qw)
54    case ('qgr', 'QGR')
55         call da_interp_lin_3d(grid%xa%qgr, iv%info(pseudo), qw)
56    end select
58    do n=iv%info(pseudo)%n1,iv%info(pseudo)%n2
59       y%pseudo(n)%u = u(1,n)
60       y%pseudo(n)%v = v(1,n)
61       y%pseudo(n)%p = p(1,n)
62       y%pseudo(n)%t = t(1,n)
63       if ( len_trim(adjustl(pseudo_var)) == 1 ) then
64          y%pseudo(n)%q = q(1,n)
65       else
66          select case(pseudo_var(1:3))
67          case ('qcw', 'QCW','qci', 'QCI','qrn', 'QRN','qsn', 'QSN', 'qgr', 'QGR')
68             y%pseudo(n)%q = qw(1,n)
69          end select
70       end if
71    end do
73    deallocate (u)
74    deallocate (v)
75    deallocate (q)
76    deallocate (p)
77    deallocate (t)
78    deallocate (qw)
80    if (trace_use_dull) call da_trace_exit("da_transform_xtoy_pseudo")
82 end subroutine da_transform_xtoy_pseudo