da_get_innov_vector_polaramv.inc
References to this file elsewhere.
1 subroutine da_get_innov_vector_polaramv( it, xb, xp, ob, iv)
2
3 !----------------------------------------------------------------------
4 ! Purpose: Calculates innovation vector, does QC for Polar AMV's
5 ! Update :
6 ! 01/24/2007 Syed RH Rizvi
7 ! Updated for "VERIFY"
8 !----------------------------------------------------------------------
9
10 implicit none
11
12 integer, intent(in) :: it ! External iteration.
13 type(xb_type), intent(in) :: xb ! first guess state.
14 type(xpose_type), intent(in) :: xp ! Domain decomposition vars.
15 type(y_type), intent(in) :: ob ! Observation structure.
16 type(ob_type), intent(inout) :: iv ! O-B structure.
17
18 integer :: n ! Loop counter.
19 integer :: i, j, k ! Index dimension.
20 integer :: num_levs ! Number of obs levels.
21
22 real :: dx, dxm ! Interpolation weights.
23 real :: dy, dym ! Interpolation weights.
24 real, dimension(1:max_ob_levels) :: model_u ! Model value u at ob location.
25 real, dimension(1:max_ob_levels) :: model_v ! Model value v at ob location.
26
27 real, dimension(xp%kms:xp%kme) :: v_p ! Model value p at ob hor. location.
28
29 integer :: itu,ituf,itvv,itvvf
30
31 if (iv % num_polaramv > 0) then
32 itu = 0; itvv = 0;
33 ituf = 0; itvvf = 0;
34
35 do n=iv%ob_numb(iv%current_ob_time-1)%polaramv + 1, iv%ob_numb(iv%current_ob_time)%polaramv
36 ! [1.3] Get horizontal interpolation weights:
37
38 num_levs = iv % polaramv(n) % info % levels
39 if (num_levs < 1) cycle
40
41 model_u(:) = 0.0
42 model_v(:) = 0.0
43
44 i = iv%polaramv(n)%loc%i
45 j = iv%polaramv(n)%loc%j
46 dx = iv%polaramv(n)%loc%dx
47 dy = iv%polaramv(n)%loc%dy
48 dxm = iv%polaramv(n)%loc%dxm
49 dym = iv%polaramv(n)%loc%dym
50
51 do k=xp%kts,xp%kte
52 v_p(k) = dym*(dxm*xb%p(i,j ,k) + dx*xb%p(i+1,j ,k)) &
53 + dy *(dxm*xb%p(i,j+1,k) + dx*xb%p(i+1,j+1,k))
54 end do
55
56 do k=1, iv % polaramv(n) % info % levels
57 iv%polaramv(n)%zk(k)=missing_r
58
59 if (iv % polaramv(n) % p(k) > 1.0) then
60 call da_to_zk(iv % polaramv(n) % p(k), v_p, xp, v_interp_p, iv%polaramv(n)%zk(k))
61 end if
62
63 if (iv%polaramv(n)%zk(k) < 0.0 .and. .not.anal_type_verify) then
64 iv % polaramv(n) % u(k) % qc = missing
65 iv % polaramv(n) % v(k) % qc = missing
66 end if
67 end do
68
69 call da_interp_lin_3d( xb % u, xp, i, j, dx, dy, dxm, dym, &
70 model_u, max_ob_levels, iv%polaramv(n)%zk, num_levs)
71 call da_interp_lin_3d( xb % v, xp, i, j, dx, dy, dxm, dym, &
72 model_v, max_ob_levels, iv%polaramv(n)%zk, num_levs)
73
74 do k = 1, iv % polaramv(n) % info % levels
75 iv % polaramv(n) % u(k) % inv = 0.0
76 iv % polaramv(n) % v(k) % inv = 0.0
77
78
79 if (ob % polaramv(n) % u(k) > missing_r .AND. &
80 iv % polaramv(n) % u(k) % qc >= obs_qc_pointer) then
81
82 iv % polaramv(n) % u(k) % inv = ob % polaramv(n) % u(k) - &
83 model_u(k)
84 end if
85
86 if (ob % polaramv(n) % v(k) > missing_r .AND. &
87 iv % polaramv(n) % v(k) % qc >= obs_qc_pointer) then
88
89 iv % polaramv(n) % v(k) % inv = ob % polaramv(n) % v(k) - &
90 model_v(k)
91 end if
92 end do
93
94 !------------------------------------------------------------------------
95 ! Perform optional maximum error check:
96 !------------------------------------------------------------------------
97
98 if (check_max_iv) call da_check_max_iv_polaramv(it, iv % polaramv(n), &
99 itu,ituf,itvv,itvvf)
100 end do
101
102 if (rootproc .and. check_max_iv_print) then
103 write(unit = check_max_iv_unit, fmt ='(A,i5,A)')'For outer iteration ',it, &
104 ', Total Rejections for Polar AMVs follows:'
105 write(unit = check_max_iv_unit, fmt = '(/,10(2(A,I6),/))') &
106 'Number of failed u-wind observations: ',ituf, ' on ',itu, &
107 'Number of failed v-wind observations: ',itvvf,' on ',itvv, &
108 'Finally Total Polar AMVs rejections ',ituf+itvvf,' on ',itu+itvv
109 end if
110 end if
111
112 end subroutine da_get_innov_vector_polaramv
113
114