da_transform_xtoy_radar.inc

References to this file elsewhere.
1 subroutine da_transform_xtoy_radar (xa, xb, iv, xp, y)
2 
3    !-----------------------------------------------------------------------
4    ! Purpose: calculate the Doppler radial velocity and 
5    ! reflectivity at the observation location from the first guess.
6    ! It is linearized. 
7    !---------------------------------------------------------------------
8  
9    implicit none
10 
11    type (x_type),  intent(in)   :: xa       ! gridded analysis increment.
12    type (xb_type), intent(in)   :: xb          ! first guess state.
13    type (ob_type), intent(in)   :: iv       ! Innovation vector (O-B).
14    type (xpose_type), intent(in):: xp       ! Domain decomposition vars.
15    type (y_type), intent(inout) :: y        ! y = h (xa) (linear)
16 
17    integer                      :: n        ! Loop counter.
18    integer                      :: i, j, k  ! Index dimension.
19    integer                      :: num_levs ! Number of obs levels.
20    real                         :: dx, dxm  ! 
21    real                         :: dy, dym  !
22 
23    real, dimension(max_ob_levels)  :: model_p  ! Model value p at ob location.
24    real, dimension(max_ob_levels)  :: model_u  ! Model value u at ob location.
25    real, dimension(max_ob_levels)  :: model_v  ! Model value v at ob location.
26    real, dimension(max_ob_levels)  :: model_w  ! Model value w at ob location.
27    real, dimension(max_ob_levels)  :: model_qrn ! Model qrn at ob location.
28    real, dimension(max_ob_levels)  :: model_qrnb! Model qrn at ob location.
29    real                            :: model_ps
30 
31    real                         :: xr,yr,zr
32 
33    if (iv%num_Radar > 0) then
34       do n=iv%ob_numb(iv%current_ob_time-1)%radar + 1, &
35          iv%ob_numb(iv%current_ob_time)%radar
36          num_levs = iv % Radar(n) % info % levels
37 
38          ! [1.3] Get dot pt. horizontal interpolation weights:
39 
40          i = iv%Radar(n)%loc%i
41          j = iv%Radar(n)%loc%j
42          dx = iv%Radar(n)%loc%dx
43          dy = iv%Radar(n)%loc%dy
44          dxm = iv%Radar(n)%loc%dxm
45          dym = iv%Radar(n)%loc%dym
46 
47          do k = 1, num_levs
48             model_qrnb(k) = iv%Radar(n)%model_qrn(k)
49             model_p(k) = iv%Radar(n)%model_p(k)
50          end do
51 
52          model_ps   = iv%Radar(n)%model_ps
53 
54          ! [1.4] Interpolate horizontally from dot points:
55          call da_interp_lin_3d(xa % u, xp, i, j, dx, dy, dxm, dym, &
56                             model_u, num_levs, iv%Radar(n)%zk, &
57                             num_levs)
58          call da_interp_lin_3d(xa % v, xp, i, j, dx, dy, dxm, dym, &
59                             model_v, num_levs, iv%Radar(n)%zk, &
60                             num_levs)
61          call da_interp_lin_3d(xa % qrn, xp, &
62                             i, j, dx, dy, dxm, dym, &
63                             model_qrn, num_levs, iv%Radar(n)%zk, &
64                             num_levs)
65          call da_interp_lin_3d(xa%wh, xp,   &
66                             i, j, dx, dy, dxm, dym, &
67                             model_w, num_levs, iv%Radar(n)%zk, &
68                             num_levs)
69 
70          ! [1.7] Calculate rv and rf at OBS location
71 
72          xr = xb%ds * (iv%Radar(n)%loc%x - iv%Radar(n)%stn_loc%x)
73          yr = xb%ds * (iv%Radar(n)%loc%y - iv%Radar(n)%stn_loc%y)
74 
75          do k = 1, num_levs
76 
77             ! y%Radar(n)%rv(k) = 0.0
78             ! y%Radar(n)%rf(k) = 0.0
79 
80             if (iv % Radar(n) % height_qc(k) /= below_model_surface .and.  &
81                  iv % Radar(n) % height_qc(k) /= above_model_lid) then
82                if (use_Radar_rv) then
83                   if (iv % Radar(n) % rv(k) % qc >= obs_qc_pointer) then
84                      zr=iv%Radar(n)%height(k) - iv%Radar(n)%stn_loc%elv
85 
86                      call da_radial_velocity_lin(y%Radar(n)%rv(k), &
87                         model_p(k), &
88                         model_u(k), model_v(k), model_w(k), model_qrn(k),    &
89                         model_ps, xr, yr, zr, model_qrnb(k))
90                   end if
91                end if
92 
93                if (use_Radar_rf) then
94                   if (iv % Radar(n) % rf(k) % qc >= obs_qc_pointer) then
95                      call da_reflectivity_lin(y%Radar(n)%rf(k), &
96                         model_qrn(k), &
97                         model_qrnb(k))
98                   end if
99                end if
100             end if
101          end do
102       end do
103    end if
104 
105 end subroutine da_transform_xtoy_radar 
106 
107