da_qc_mhs.inc
References to this file elsewhere.
1 subroutine da_qc_mhs (i, nchan, ob, iv)
2
3 !---------------------------------------------------------------------------
4 ! Purpose: perform quality control for mhs data.
5 !---------------------------------------------------------------------------
6
7 implicit none
8
9 integer, intent(in) :: i ! sensor index.
10 integer, intent(in) :: nchan ! number of channel
11 type (y_type), intent(in) :: ob ! Observation structure.
12 type (iv_type), intent(inout) :: iv ! O-B structure.
13
14 ! local variables
15 integer :: n,scanpos,k,isflg,ios,fgat_rad_unit
16 logical :: lmix
17 real :: si
18 ! real :: satzen
19 integer :: ngood(nchan),nrej(nchan),nrej_omb_abs(nchan), &
20 nrej_omb_std(nchan), &
21 nrej_mixsurface,nrej_windowchanl, nrej_si, &
22 nrej_clw,nrej_topo, num_proc_domain, &
23 nrej_limb
24
25 character(len=30) :: filename
26
27 if (trace_use) call da_trace_entry("da_qc_mhs")
28
29 ngood(:) = 0
30 nrej(:) = 0
31 nrej_omb_abs(:) = 0
32 nrej_omb_std(:) = 0
33 nrej_mixsurface = 0
34 nrej_windowchanl= 0
35 nrej_si = 0
36 nrej_clw = 0
37 nrej_topo = 0
38 nrej_limb = 0
39
40 num_proc_domain = 0
41
42 do n= iv%instid(i)%info%n1,iv%instid(i)%info%n2
43
44 if (iv%instid(i)%info%proc_domain(1,n)) &
45 num_proc_domain = num_proc_domain + 1
46
47 ! 0.0 initialise QC flags by assuming good obs
48 !---------------------------------------------
49 iv%instid(i)%tb_qc(:,n) = qc_good
50
51 ! a. reject all channels over mixture surface type
52 !------------------------------------------------------
53 isflg = iv%instid(i)%isflg(n)
54 lmix = (isflg==4) .or. (isflg==5) .or. (isflg==6) .or. (isflg==7)
55 if (lmix) then
56 iv%instid(i)%tb_qc(:,n) = qc_bad
57 if (iv%instid(i)%info%proc_domain(1,n)) &
58 nrej_mixsurface = nrej_mixsurface + 1
59 end if
60
61 ! b. reject channels 1~2 over land/sea-ice/snow
62 !------------------------------------------------------
63 if (isflg > 0) then
64 iv%instid(i)%tb_qc(1:2,n) = qc_bad
65 if (iv%instid(i)%info%proc_domain(1,n)) &
66 nrej_windowchanl = nrej_windowchanl + 1
67 end if
68
69 ! reject limb obs
70 !------------------------------------------------------
71 scanpos = iv%instid(i)%scanpos(n)
72 if (scanpos <= 8 .or. scanpos >= 83) then
73 iv%instid(i)%tb_qc(:,n) = qc_bad
74 if (iv%instid(i)%info%proc_domain(1,n)) &
75 nrej_limb = nrej_limb + 1
76 end if
77
78 ! satzen = rad%satzen
79 ! if (abs(satzen) > 45.0) &
80 ! iv%instid(i)%tb_qc(:,n) = qc_bad
81
82 ! d. check cloud/precipitation
83 !-----------------------------------------------------------
84 if (ob%instid(i)%tb(1,n) > 0.0 .and. &
85 ob%instid(i)%tb(2,n) > 0.0) then
86 si = ob%instid(i)%tb(1,n) - ob%instid(i)%tb(2,n)
87 if (si >= 3.0) then
88 iv%instid(i)%tb_qc(:,n) = qc_bad
89 iv%instid(i)%cloud_flag(:,n) = qc_bad
90 if (iv%instid(i)%info%proc_domain(1,n)) &
91 nrej_si = nrej_si + 1
92 end if
93 end if
94
95 if (iv%instid(i)%clwp(n) >= 0.2) then
96 iv%instid(i)%tb_qc(:,n) = qc_bad
97 iv%instid(i)%cloud_flag(:,n) = qc_bad
98 if (iv%instid(i)%info%proc_domain(1,n)) &
99 nrej_clw = nrej_clw + 1
100 end if
101
102 ! e. check surface height/pressure
103 !------------------------------------------------------
104 ! sfchgt = ivrad%info(n)%elv
105 ! if (sfchgt >=) then
106 !
107 ! else
108 ! end if
109
110 if ((isflg .ne. 0) .and. (iv%instid(i)%ps(n) < 800.0)) then
111 iv%instid(i)%tb_qc(5,n) = qc_bad
112 if (iv%instid(i)%info%proc_domain(1,n)) &
113 nrej_topo = nrej_topo + 1
114 end if
115
116 ! g. check iuse (pre-rejected channels by .info files)
117 !------------------------------------------------------
118 do k = 1, nchan
119 if (satinfo(i)%iuse(k) .eq. -1) &
120 iv%instid(i)%tb_qc(k,n) = qc_bad
121 end do
122
123 ! f. check innovation
124 !-----------------------------------------------------------
125 do k = 1, nchan
126 ! absolute departure check
127 if (abs(iv%instid(i)%tb_inv(k,n)) > 15.0) then
128 iv%instid(i)%tb_qc(k,n) = qc_bad
129 if (iv%instid(i)%info%proc_domain(1,n)) &
130 nrej_omb_abs(k) = nrej_omb_abs(k) + 1
131 end if
132
133 ! relative departure check
134 if (use_error_factor_rad) then
135 iv%instid(i)%tb_error(k,n) = &
136 satinfo(i)%error_std(k)*satinfo(i)%error_factor(k)
137 else
138 iv%instid(i)%tb_error(k,n) = satinfo(i)%error_std(k)
139 end if
140
141 if (abs(iv%instid(i)%tb_inv(k,n)) > 3.0*iv%instid(i)%tb_error(k,n)) then
142 iv%instid(i)%tb_qc(k,n) = qc_bad
143 if (iv%instid(i)%info%proc_domain(1,n)) &
144 nrej_omb_std(k) = nrej_omb_std(k) + 1
145 end if
146
147 ! final QC decision
148 if (iv%instid(i)%tb_qc(k,n) == qc_bad) then
149 iv%instid(i)%tb_error(k,n) = 500.0
150 if (iv%instid(i)%info%proc_domain(1,n)) &
151 nrej(k) = nrej(k) + 1
152 else
153 if (iv%instid(i)%info%proc_domain(1,n)) &
154 ngood(k) = ngood(k) + 1
155 end if
156 end do ! chan
157 end do ! end loop pixel
158
159 ! Do inter-processor communication to gather statistics.
160
161 call da_proc_sum_int (num_proc_domain)
162 call da_proc_sum_int (nrej_mixsurface)
163 call da_proc_sum_int (nrej_windowchanl)
164 call da_proc_sum_int (nrej_si)
165 call da_proc_sum_int (nrej_clw)
166 call da_proc_sum_int (nrej_topo)
167 call da_proc_sum_int (nrej_limb)
168 call da_proc_sum_ints (nrej_omb_abs(:))
169 call da_proc_sum_ints (nrej_omb_std(:))
170 call da_proc_sum_ints (nrej(:))
171 call da_proc_sum_ints (ngood(:))
172
173 if (rootproc) then
174 if (num_fgat_time > 1) then
175 write(filename,'(a,i2.2)') 'qcstat_'//trim(iv%instid(i)%rttovid_string)//'_',iv%time
176 else
177 filename = 'qcstat_'//trim(iv%instid(i)%rttovid_string)
178 end if
179 call da_get_unit(fgat_rad_unit)
180 open(fgat_rad_unit,file=trim(filename),form='formatted',iostat=ios)
181 if (ios /= 0) then
182 write(unit=message(1),fmt='(A,A)') 'error opening the output file ', filename
183 call da_error(__FILE__,__LINE__,message(1:1))
184 end if
185
186 write(fgat_rad_unit, fmt='(/a/)') &
187 'Quality Control Statistics for '//iv%instid(i)%rttovid_string
188 write(fgat_rad_unit,'(a20,i7)') ' num_proc_domain = ', num_proc_domain
189 write(fgat_rad_unit,'(a20,i7)') ' nrej_mixsurface = ', nrej_mixsurface
190 write(fgat_rad_unit,'(a20,i7)') ' nrej_windowchanl = ', nrej_windowchanl
191 write(fgat_rad_unit,'(a20,i7)') ' nrej_si = ', nrej_si
192 write(fgat_rad_unit,'(a20,i7)') ' nrej_clw = ', nrej_clw
193 write(fgat_rad_unit,'(a20,i7)') ' nrej_topo = ', nrej_topo
194 write(fgat_rad_unit,'(a20,i7)') ' nrej_limb = ', nrej_limb
195 write(fgat_rad_unit,'(a20)') ' nrej_omb_abs(:) = '
196 write(fgat_rad_unit,'(10i7)') nrej_omb_abs(:)
197 write(fgat_rad_unit,'(a20)') ' nrej_omb_std(:) = '
198 write(fgat_rad_unit,'(10i7)') nrej_omb_std(:)
199 write(fgat_rad_unit,'(a20)') ' nrej(:) = '
200 write(fgat_rad_unit,'(10i7)') nrej(:)
201 write(fgat_rad_unit,'(a20)') ' ngood(:) = '
202 write(fgat_rad_unit,'(10i7)') ngood(:)
203
204 close(fgat_rad_unit)
205 call da_free_unit(fgat_rad_unit)
206 end if
207
208 if (trace_use) call da_trace_exit("da_qc_mhs")
209
210 end subroutine da_qc_mhs
211
212