1 subroutine da_get_innov_vector_ssmi_rv (it,num_qcstat_conv, grid, ob, iv)
3 !-----------------------------------------------------------------------
5 ! Updated for Analysis on Arakawa-C grid
6 ! Author: Syed RH Rizvi, MMM/ESSL/NCAR, Date: 10/22/2008
7 !-----------------------------------------------------------------------
11 integer, intent(in) :: it ! External iteration.
12 type(domain), intent(in) :: grid ! first guess state.
13 type(y_type), intent(in) :: ob ! Observation structure.
14 type(iv_type), intent(inout) :: iv ! O-B structure.
15 integer, intent(inout) :: num_qcstat_conv(:,:,:,:)
17 integer :: n ! Loop counter.
18 integer :: i, j ! Index dimension.
19 real :: dx, dxm ! Interpolation weights.
20 real :: dy, dym ! Interpolation weights.
21 real :: model_tpw ! Model value tpw at oblocation
22 real :: model_speed ! Model value speed at oblocation
23 integer :: itpw, itpwf, ispeed, ispeedf
25 if (trace_use) call da_trace_entry("da_get_innov_vector_ssmi_rv")
28 do n=iv%info(ssmi_rv)%n1,iv%info(ssmi_rv)%n2
29 if ( iv % ssmi_rv(n) % speed % qc == fails_error_max) iv % ssmi_rv(n) % speed % qc = 0
30 if ( iv % ssmi_rv(n) % tpw % qc == fails_error_max) iv % ssmi_rv(n) % tpw % qc = 0
34 do n=iv%info(ssmi_rv)%n1,iv%info(ssmi_rv)%n2
36 ! compute innovation vector
37 ! =========================
39 ! Obs coordinates on model grid
43 i = iv%info(ssmi_rv)%i(1,n)
44 j = iv%info(ssmi_rv)%j(1,n)
45 dx = iv%info(ssmi_rv)%dx(1,n)
46 dy = iv%info(ssmi_rv)%dy(1,n)
47 dxm = iv%info(ssmi_rv)%dxm(1,n)
48 dym = iv%info(ssmi_rv)%dym(1,n)
50 iv % ssmi_rv(n) % tpw % inv = 0.0
51 if (abs(ob % ssmi_rv(n) % tpw - missing_r) > 1.0 .and. iv % ssmi_rv(n) % tpw % qc >= obs_qc_pointer) then
52 model_tpw = dym*(dxm*grid%xb%tpw(i,j) + dx*grid%xb%tpw(i+1,j)) + dy *(dxm*grid%xb%tpw(i,j+1) + dx*grid%xb%tpw(i+1,j+1))
53 iv % ssmi_rv(n) % tpw % inv = ob % ssmi_rv(n) % tpw - model_tpw
58 iv % ssmi_rv(n) % speed % inv = 0.0
59 if (abs(ob % ssmi_rv(n) % speed - missing_r) > 1.0 .and. &
60 iv % ssmi_rv(n) % speed % qc >= obs_qc_pointer) then
62 model_speed = dym*(dxm*grid%xb%speed(i,j ) + dx*grid%xb%speed(i+1,j )) &
63 + dy *(dxm*grid%xb%speed(i,j+1) + dx*grid%xb%speed(i+1,j+1))
64 iv % ssmi_rv(n) % speed % inv = ob % ssmi_rv(n) % speed - model_speed
68 !------------------------------------------------------------------
69 ! Perform optional maximum error check:
70 !------------------------------------------------------------------
73 call da_check_max_iv_ssmi_rv(iv, it, num_qcstat_conv)
75 if (trace_use) call da_trace_exit("da_get_innov_vector_ssmi_rv")
77 end subroutine da_get_innov_vector_ssmi_rv