da_qc_airs.inc
References to this file elsewhere.
1 subroutine da_qc_airs (i, nchan, ob, iv)
2
3 !---------------------------------------------------------------------------
4 ! Purpose: perform quality control for AQUA/EOS-2-AIRS 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 (ob_type), intent(inout) :: iv ! O-B structure.
13
14
15 ! local variables
16 integer :: n,k,isflg,ios,fgat_rad_unit
17 ! integer :: scanpos
18 ! logical :: lmix
19 ! real :: satzen
20 !!
21 !! TV change: added new diagnostic local varaibles
22 integer :: ngood(nchan),nrej(nchan),nrej_omb_abs(nchan), &
23 nrej_omb_std(nchan),nrej_limb, &
24 nrej_landsurface,nrej_windowchshort,nrej_windowchlong, &
25 nrej_clw,nrej_sst,nrej_topo,num_rad, num_proc_domain
26 ! integer :: nrej_si
27 !
28 ! integer :: ngood(nchan),nrej(nchan),nrej_omb_abs(nchan), &
29 ! nrej_omb_std(nchan), &
30 ! nrej_mixsurface,nrej_windowchanl, nrej_si, &
31 ! nrej_clw,nrej_topo,num_rad, num_proc_domain, &
32 ! nrej_limb
33 !
34 real :: SST_model, SST_airs, diffSST
35
36 !! TV end change
37
38
39 character(len=30) :: filename
40
41 if (trace_use) call da_trace_entry("da_qc_airs")
42
43 !! TV change: new diagnosic locals
44 ngood(:) = 0
45 nrej(:) = 0
46 nrej_omb_abs(:) = 0
47 nrej_omb_std(:) = 0
48 nrej_landsurface = 0
49 nrej_windowchshort= 0
50 nrej_windowchlong= 0
51 nrej_sst= 0
52 ! nrej_si = 0
53 nrej_clw = 0
54 nrej_topo = 0
55
56 ! nrej_mixsurface = 0
57 ! nrej_windowchanl= 0
58
59 !! TV end change
60
61 nrej_limb = 0
62 num_rad = iv%ob_numb(iv%current_ob_time)%radiance(i)- &
63 iv%ob_numb(iv%current_ob_time-1)%radiance(i)
64 ! num_rad = iv%instid(i)%num_rad
65 num_proc_domain = 0
66
67 if (num_rad > 0) then
68
69 ! do n= 1, iv%instid(i)%num_rad ! loop for pixel
70 do n= iv%ob_numb(iv%current_ob_time-1)%radiance(i)+1, iv%ob_numb(iv%current_ob_time)%radiance(i)
71
72 if (iv%instid(i)%proc_domain(n)) &
73 num_proc_domain = num_proc_domain + 1
74
75 ! 0.0 initialise QC by flags assuming good obs
76 !---------------------------------------------
77 iv%instid(i)%tb_qc(:,n) = qc_good
78
79 !! TV change: include sea locations only
80 ! a. reject all channels over land and mixture
81 !------------------------------------------------------
82 isflg = iv%instid(i)%isflg(n)
83 if (isflg > 0) then
84 iv%instid(i)%tb_qc(:,n) = -1
85 if (iv%instid(i)%proc_domain(n)) &
86 nrej_landsurface = nrej_landsurface + 1
87 end if
88
89 ! ! a. reject all channels over mixture surface type
90 ! !------------------------------------------------------
91 ! isflg = iv%instid(i)%isflg(n)
92 ! lmix = (isflg==4) .or. (isflg==5) .or. (isflg==6) .or. (isflg==7)
93 ! if (lmix) then
94 ! iv%instid(i)%tb_qc(:,n) = qc_bad
95 ! if (iv%instid(i)%proc_domain(n)) &
96 ! nrej_mixsurface = nrej_mixsurface + 1
97 ! end if
98
99 ! ! b. reject all channels over land/sea-ice/snow
100 ! !------------------------------------------------------
101 ! if (isflg > 0) then
102 ! iv%instid(i)%tb_qc(:,n) = qc_bad
103 ! if (iv%instid(i)%proc_domain(n)) &
104 ! nrej_windowchanl = nrej_windowchanl + 1
105 ! end if
106
107 !! TV end change: include sea locations only
108
109
110 ! reject limb obs
111 !------------------------------------------------------
112 !scanpos = iv%instid(i)%scanpos(n)
113 !if (scanpos <= 3 .or. scanpos >= 88) then
114 ! iv%instid(i)%tb_qc(:,n) = qc_bad
115 ! if (iv%instid(i)%proc_domain(n)) &
116 ! nrej_limb = nrej_limb + 1
117 !end if
118
119
120 !! TV change: QC for clouds
121
122 ! c. Check for model clouds!!!!!!
123 !-----------------------------------------------------------
124
125 if (iv%instid(i)%clwp(n) > 0.05) then
126 iv%instid(i)%tb_qc(:,n) = qc_bad
127 iv%instid(i)%cloud_flag(:,n) = -1
128 if (iv%instid(i)%proc_domain(n)) &
129 nrej_clw = nrej_clw + 1
130 end if
131
132 ! d. Crude check for clouds in obs (assuming obs are used over ocean only)
133 ! Use long wave window channel #914 - 10.662 nm (965.43 cm^-1)
134 ! should be warmer than freezing temperature of the sea
135 !-----------------------------------------------------------
136 !
137 if(ob%instid(i)%tb(129,n) < 271.) then
138 iv%instid(i)%tb_qc(:,n) = qc_bad
139 if (iv%instid(i)%proc_domain(n)) &
140 nrej_windowchlong = nrej_windowchlong + 1
141 end if
142
143 ! e. Check for contaminated obs in warmest near-infrared: Sun contamination during day
144 !-----------------------------------------------------------
145 !
146 SST_airs=ob%instid(i)%tb(272,n) !! short wave window channel #2328 - 3.882 nm (2616.38 cm^-1)
147 if(SST_airs > 307.) then
148 iv%instid(i)%tb_qc(257:281,n) = qc_bad
149 if (iv%instid(i)%proc_domain(n)) &
150 nrej_windowchshort = nrej_windowchshort + 1
151 end if
152
153
154 ! f. Check for cloud free in obs (assuming obs are used over ocean only)
155 ! Criterion: model SST within 2 K of transparent (hottest) short wave window channel
156 ! includes check for contaminated near-infrared
157 !-----------------------------------------------------------
158 !
159 SST_model=iv%instid(i)%ts(n) !! SST in the model
160 diffSST=abs(SST_model-SST_airs)
161 if(SST_airs < 307 .and. diffSST > 2.) then
162 iv%instid(i)%tb_qc(:,n) = qc_bad
163 if (iv%instid(i)%proc_domain(n)) &
164 nrej_sst = nrej_sst + 1
165 end if
166
167 !! TV change
168
169 ! e. check surface height/pressure
170 !-----------------------------------------------------------
171 ! sfchgt = ivrad%info(n)%elv
172 ! if (sfchgt >=) then
173 ! else
174 ! end if
175
176 !if ((isflg .ne. 0) .and. (iv%instid(i)%ps(n) < 850.)) then
177 ! iv%instid(i)%tb_qc(5,n) = qc_bad
178 ! if (iv%instid(i)%proc_domain(n)) &
179 ! nrej_topo = nrej_topo + 1
180 !end if
181
182 ! g. check iuse from information file (channel selection)
183 !-----------------------------------------------------------
184 do k = 1, nchan
185 if (satinfo(i)%iuse(k) .eq. -1) &
186 iv%instid(i)%tb_qc(k,n) = qc_bad
187 end do
188
189 ! f. check innovation
190 !-----------------------------------------------------------
191 do k = 1, nchan
192
193 ! 1. check absolute value of innovation
194 !------------------------------------------------
195 if (abs(iv%instid(i)%tb_inv(k,n)) > 15.) then
196 iv%instid(i)%tb_qc(k,n) = qc_bad
197 if (iv%instid(i)%proc_domain(n)) &
198 nrej_omb_abs(k) = nrej_omb_abs(k) + 1
199 end if
200
201 ! 2. check relative value of innovation
202 ! and assign of the observation error (standard deviation)
203 !------------------------------------------------------------------------
204 if (use_error_factor_rad) then ! if use error tuning factor
205 iv%instid(i)%tb_error(k,n) = &
206 satinfo(i)%error(k)*satinfo(i)%error_factor(k)
207 else
208 iv%instid(i)%tb_error(k,n) = satinfo(i)%error(k)
209 end if
210
211 if (abs(iv%instid(i)%tb_inv(k,n)) > 3.0*iv%instid(i)%tb_error(k,n)) then
212 iv%instid(i)%tb_qc(k,n) = qc_bad
213 if (iv%instid(i)%proc_domain(n)) &
214 nrej_omb_std(k) = nrej_omb_std(k) + 1
215 end if
216
217 ! 3. Final QC decision
218 !---------------------------------------------
219 if (iv%instid(i)%tb_qc(k,n) == qc_bad) then ! bad obs
220 iv%instid(i)%tb_error(k,n) = 500.0
221 if (iv%instid(i)%proc_domain(n)) &
222 nrej(k) = nrej(k) + 1
223 else ! good obs
224 if (iv%instid(i)%proc_domain(n)) &
225 ngood(k) = ngood(k) + 1
226 end if
227 end do ! chan
228 end do ! end loop pixel
229 end if
230
231 ! Do inter-processor communication to gather statistics.
232 call da_proc_sum_int (num_proc_domain)
233 !! TV change
234 call da_proc_sum_int ( nrej_landsurface )
235 call da_proc_sum_int ( nrej_windowchlong)
236 call da_proc_sum_int ( nrej_windowchshort)
237 call da_proc_sum_int ( nrej_sst)
238 ! call da_proc_sum_int ( nrej_si )
239 !! TV end change
240 call da_proc_sum_int ( nrej_clw )
241 call da_proc_sum_int ( nrej_topo )
242 call da_proc_sum_int (nrej_limb)
243 call da_proc_sum_ints (nrej_omb_abs(:))
244 call da_proc_sum_ints (nrej_omb_std(:))
245 call da_proc_sum_ints (nrej(:))
246 call da_proc_sum_ints (ngood(:))
247
248 if (rootproc) then
249 if (num_fgat_time > 1) then
250 write(filename,'(a,i2.2)') 'qcstat_'//trim(iv%instid(i)%rttovid_string)//'_',iv%current_ob_time
251 else
252 filename = 'qcstat_'//trim(iv%instid(i)%rttovid_string)
253 end if
254
255 call da_get_unit(fgat_rad_unit)
256 open(fgat_rad_unit,file=trim(filename),form='formatted',iostat=ios)
257 if (ios /= 0) Then
258 write(unit=message(1),fmt='(A,A)') 'error opening the output file ', filename
259 call da_error(__FILE__,__LINE__,message(1:1))
260 end if
261
262 write(fgat_rad_unit, fmt='(/a/)') ' Quality Control Statistics for '//iv%instid(i)%rttovid_string
263 write(fgat_rad_unit,'(a20,i7)') ' num_proc_domain = ', num_proc_domain
264 !!TV change
265 write(fgat_rad_unit,'(a20,i7)') ' nrej_landsurface = ', nrej_landsurface
266 write(fgat_rad_unit,'(a20,i7)') ' nrej_windowchlong = ', nrej_windowchlong
267 write(fgat_rad_unit,'(a20,i7)') ' nrej_windowchshort = ', nrej_windowchshort
268 write(fgat_rad_unit,'(a20,i7)') ' nrej_sst = ', nrej_sst
269 ! write(fgat_rad_unit,'(a20,i7)') ' nrej_mixsurface = ', nrej_mixsurface
270 ! write(fgat_rad_unit,'(a20,i7)') ' nrej_windowchanl = ', nrej_windowchanl
271 ! write(fgat_rad_unit,'(a20,i7)') ' nrej_si = ', nrej_si
272 !!TV end change
273 write(fgat_rad_unit,'(a20,i7)') ' nrej_clw = ', nrej_clw
274 write(fgat_rad_unit,'(a20,i7)') ' nrej_topo = ', nrej_topo
275 write(fgat_rad_unit,'(a20,i7)') ' nrej_limb = ', nrej_limb
276 write(fgat_rad_unit,'(a20)') ' nrej_omb_abs(:) = '
277 write(fgat_rad_unit,'(10i7)') nrej_omb_abs(:)
278 write(fgat_rad_unit,'(a20)') ' nrej_omb_std(:) = '
279 write(fgat_rad_unit,'(10i7)') nrej_omb_std(:)
280 write(fgat_rad_unit,'(a20)') ' nrej(:) = '
281 write(fgat_rad_unit,'(10i7)') nrej(:)
282 write(fgat_rad_unit,'(a20)') ' ngood(:) = '
283 write(fgat_rad_unit,'(10i7)') ngood(:)
284
285 close(fgat_rad_unit)
286 call da_free_unit(fgat_rad_unit)
287 end if
288
289 if (trace_use) call da_trace_exit("da_qc_airs")
290
291 end subroutine da_qc_airs
292