da_transform_xtoy_radar.inc

References to this file elsewhere.
1 subroutine da_transform_xtoy_radar (grid, iv, 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 (domain),  intent(in)    :: grid
12    type (iv_type), intent(in)    :: iv       ! Innovation vector (O-B).
13    type (y_type),  intent(inout) :: y        ! y = h (grid%xa) (linear)
14 
15    integer :: n                   ! Loop counter.      
16    integer :: k                         ! Index dimension.   
17 
18    real, allocatable :: model_p(:,:)
19    real, allocatable :: model_u(:,:)
20    real, allocatable :: model_v(:,:)
21    real, allocatable :: model_w(:,:)
22    real, allocatable :: model_qrn(:,:)
23    real, allocatable :: model_qrnb(:,:)
24    real, allocatable :: model_ps(:)
25    real    :: xr,yr,zr
26 
27    real    :: alog_10
28 
29    if (trace_use) call da_trace_entry("da_transform_xtoy_radar")
30 
31    alog_10 = alog(10.0)
32 
33    allocate (model_p(iv%info(radar)%max_lev,iv%info(radar)%n1:iv%info(radar)%n2))
34    allocate (model_u(iv%info(radar)%max_lev,iv%info(radar)%n1:iv%info(radar)%n2))
35    allocate (model_v(iv%info(radar)%max_lev,iv%info(radar)%n1:iv%info(radar)%n2))
36    allocate (model_w(iv%info(radar)%max_lev,iv%info(radar)%n1:iv%info(radar)%n2))
37    allocate (model_qrn(iv%info(radar)%max_lev,iv%info(radar)%n1:iv%info(radar)%n2))
38    allocate (model_qrnb(iv%info(radar)%max_lev,iv%info(radar)%n1:iv%info(radar)%n2))
39    allocate (model_ps(iv%info(radar)%n1:iv%info(radar)%n2))
40 
41    do n=iv%info(radar)%n1,iv%info(radar)%n2
42       do k = 1, iv%info(radar)%levels(n)
43          model_qrnb(k,n) = iv%radar(n)%model_qrn(k)
44          model_p(k,n)    = iv%radar(n)%model_p(k)
45       end do
46 
47       model_ps(n) = iv%radar(n)%model_ps
48    end do
49 
50    ! [1.4] Interpolate horizontally from dot points:
51    call da_interp_lin_3d (grid%xa%u,   iv%info(radar), model_u)
52    call da_interp_lin_3d (grid%xa%v,   iv%info(radar), model_v)
53    call da_interp_lin_3d (grid%xa%qrn, iv%info(radar), model_qrn)
54    call da_interp_lin_3d (grid%xa%wh,  iv%info(radar), model_w)
55 
56    do n=iv%info(radar)%n1,iv%info(radar)%n2
57 
58       ! [1.7] Calculate rv and rf at OBS location
59 
60       xr = grid%xb%ds * (iv%info(radar)%x(1,n) - iv%radar(n)%stn_loc%x)
61       yr = grid%xb%ds * (iv%info(radar)%y(1,n) - iv%radar(n)%stn_loc%y)
62 
63       do k = 1, iv%info(radar)%levels(n)
64          if (iv % radar(n) % height_qc(k) /= below_model_surface .and.  &
65               iv % radar(n) % height_qc(k) /= above_model_lid) then
66             if (use_radar_rv) then
67                if (iv % radar(n) % rv(k) % qc >= obs_qc_pointer) then
68                   zr=iv%radar(n)%height(k) - iv%radar(n)%stn_loc%elv
69 
70                   call da_radial_velocity_lin(y%radar(n)%rv(k), &
71                      model_p(k,n), &
72                      model_u(k,n), model_v(k,n), model_w(k,n), model_qrn(k,n),    &
73                      model_ps(n), xr, yr, zr, model_qrnb(k,n))
74                end if
75             end if
76 
77             if (use_radar_rf) then
78                if (iv % radar(n) % rf(k) % qc >= obs_qc_pointer) then
79                   y%radar(n)%rf(k) = leh2 * model_qrn(k,n) /(model_qrnb(k,n)*alog_10) 
80                end if
81             end if
82          end if
83       end do
84    end do
85 
86    deallocate (model_p)
87    deallocate (model_u)
88    deallocate (model_v)
89    deallocate (model_w)
90    deallocate (model_qrn)
91    deallocate (model_qrnb)
92    deallocate (model_ps)
93 
94    if (trace_use) call da_trace_exit("da_transform_xtoy_radar")
95 
96 end subroutine da_transform_xtoy_radar 
97 
98