1 subroutine da_qc_mwts (it, i, nchan, ob, iv)
3 !---------------------------------------------------------------------------
4 ! Purpose: perform quality control for fy3 mwts data.
5 ! Dong Peiming modified from da_qc_amsua 2012/03/15
6 !---------------------------------------------------------------------------
10 integer, intent(in) :: it ! outer loop count
11 integer, intent(in) :: i ! sensor index.
12 integer, intent(in) :: nchan ! number of channel
13 type (y_type), intent(in) :: ob ! Observation structure.
14 type (iv_type), intent(inout) :: iv ! O-B structure.
18 integer :: n,scanpos,k,isflg,ios,fgat_rad_unit
22 integer :: ngood(nchan),nrej(nchan),nrej_omb_abs(nchan), &
23 nrej_omb_std(nchan), &
24 nrej_mixsurface,nrej_windowchanl, nrej_si, &
25 nrej_clw,nrej_topo, num_proc_domain, &
28 character(len=30) :: filename
30 if (trace_use) call da_trace_entry("da_qc_mwts")
45 do n= iv%instid(i)%info%n1,iv%instid(i)%info%n2
47 if (iv%instid(i)%info%proc_domain(1,n)) &
48 num_proc_domain = num_proc_domain + 1
50 ! 0.0 initialise QC by flags assuming good obs
51 !---------------------------------------------
52 iv%instid(i)%tb_qc(:,n) = qc_good
54 ! a. reject all channels over mixture surface type
55 !------------------------------------------------------
56 isflg = iv%instid(i)%isflg(n)
57 lmix = (isflg==4) .or. (isflg==5) .or. (isflg==6) .or. (isflg==7)
59 iv%instid(i)%tb_qc(:,n) = qc_bad
60 if (iv%instid(i)%info%proc_domain(1,n)) &
61 nrej_mixsurface = nrej_mixsurface + 1
63 ! b. reject channels 1 over land/sea-ice/snow
64 !------------------------------------------------------
66 iv%instid(i)%tb_qc(1:1,n) = qc_bad
67 if (iv%instid(i)%info%proc_domain(1,n)) &
68 nrej_windowchanl = nrej_windowchanl + 1
69 ! reject whole pixel if not over sea for global case
70 if (global) iv%instid(i)%tb_qc(:,n) = qc_bad
71 if (only_sea_rad) iv%instid(i)%tb_qc(:,n) = qc_bad
74 ! c. reject channels 13,14(above top model 10mb),15
75 !------------------------------------------------------
76 !Dongpm iv%instid(i)%tb_qc(13:15,n) = qc_bad
79 !------------------------------------------------------
80 scanpos = iv%instid(i)%scanpos(n)
81 if (scanpos <= 2 .or. scanpos >= 14) then
82 iv%instid(i)%tb_qc(:,n) = qc_bad
83 if (iv%instid(i)%info%proc_domain(1,n)) &
84 nrej_limb = nrej_limb + 1
88 ! if (abs(satzen) > 45.0) iv%instid(i)%tb_qc(:,n) = qc_bad
90 ! d. check precipitation
91 !-----------------------------------------------------------
92 !Dongpm if (ob%instid(i)%tb(1,n) > 0.0 .and. &
93 !Dongpm ob%instid(i)%tb(15,n) > 0.0) then
94 !Dongpm si = ob%instid(i)%tb(1,n) - ob%instid(i)%tb(15,n)
95 !Dongpm if (si >= 3.0) then
96 !Dongpm iv%instid(i)%tb_qc(:,n) = qc_bad
97 !Dongpm iv%instid(i)%cloud_flag(:,n) = qc_bad
98 !Dongpm if (iv%instid(i)%info%proc_domain(1,n)) &
99 !Dongpm nrej_si = nrej_si + 1
103 if (ob%instid(i)%tb(1,n) > 0.0) then
104 si = iv%instid(i)%tb_inv(1,n)
105 if (isflg .eq.0 .AND. si >= 3.0) then
106 iv%instid(i)%tb_qc(1:3,n) = qc_bad
107 iv%instid(i)%cloud_flag(1:3,n) = qc_bad
108 if (iv%instid(i)%info%proc_domain(1,n)) &
109 nrej_si = nrej_si + 1
110 elseif (isflg .gt.0 .AND. si >= 1.5) then
111 iv%instid(i)%tb_qc(1:3,n) = qc_bad
112 iv%instid(i)%cloud_flag(1:3,n) = qc_bad
113 if (iv%instid(i)%info%proc_domain(1,n)) &
114 nrej_si = nrej_si + 1
118 if (iv%instid(i)%clwp(n) >= 0.2) then
119 iv%instid(i)%tb_qc(:,n) = qc_bad
120 iv%instid(i)%cloud_flag(:,n) = qc_bad
121 if (iv%instid(i)%info%proc_domain(1,n)) &
122 nrej_clw = nrej_clw + 1
125 ! 3.1 Estimate Cloud Liquid Water (CLW) in mm over sea
126 ! (Grody etal. 2001, JGR, Equation 5b,7c,7d,9)
127 !---------------------------------------------------------
128 ! if (isflg == 0) then
129 ! coszen = cos(iv%instid(i)%satzen(n))
130 ! d0 = 8.24-(2.622-1.846*coszen)*coszen
133 ! ts = iv%instid(i)%ts(n)
134 ! tb1 = ob%instid(i)%tb(1,n)
135 ! tb2 = ob%instid(i)%tb(2,n)
136 ! clw = coszen*(d0+d1*log(ts-tb1)+d2*log(ts-tb2))
141 ! e. check surface height/pressure
142 !-----------------------------------------------------------
143 ! sfchgt = ivrad%info(n)%elv
144 ! if (sfchgt >=) then
148 if ((isflg .ne. 0) .and. (iv%instid(i)%ps(n) < 850.0)) then
149 iv%instid(i)%tb_qc(2,n) = qc_bad
150 if (iv%instid(i)%info%proc_domain(1,n)) &
151 nrej_topo = nrej_topo + 1
155 !-----------------------------------------------------------
157 if (satinfo(i)%iuse(k) .eq. -1) &
158 iv%instid(i)%tb_qc(k,n) = qc_bad
161 ! f. check innovation
162 !-----------------------------------------------------------
165 ! absolute departure check
166 if (abs(iv%instid(i)%tb_inv(k,n)) > 15.0) then
167 iv%instid(i)%tb_qc(k,n) = qc_bad
168 if (iv%instid(i)%info%proc_domain(1,n)) &
169 nrej_omb_abs(k) = nrej_omb_abs(k) + 1
172 ! relative departure check
173 if (use_error_factor_rad) then
174 iv%instid(i)%tb_error(k,n) = &
175 satinfo(i)%error_std(k)*satinfo(i)%error_factor(k)
177 iv%instid(i)%tb_error(k,n) = satinfo(i)%error_std(k)
180 if (abs(iv%instid(i)%tb_inv(k,n)) > 3.0*iv%instid(i)%tb_error(k,n)) then
181 iv%instid(i)%tb_qc(k,n) = qc_bad
182 if (iv%instid(i)%info%proc_domain(1,n)) &
183 nrej_omb_std(k) = nrej_omb_std(k) + 1
187 if (iv%instid(i)%tb_qc(k,n) == qc_bad) then
188 iv%instid(i)%tb_error(k,n) = 500.0
189 if (iv%instid(i)%info%proc_domain(1,n)) &
190 nrej(k) = nrej(k) + 1
192 if (iv%instid(i)%info%proc_domain(1,n)) &
193 ngood(k) = ngood(k) + 1
197 end do ! end loop pixel
199 ! Do inter-processor communication to gather statistics.
200 call da_proc_sum_int (num_proc_domain)
201 call da_proc_sum_int (nrej_mixsurface)
202 call da_proc_sum_int (nrej_windowchanl)
203 call da_proc_sum_int (nrej_si )
204 call da_proc_sum_int (nrej_clw)
205 call da_proc_sum_int (nrej_topo)
206 call da_proc_sum_int (nrej_limb)
207 call da_proc_sum_ints (nrej_omb_abs(:))
208 call da_proc_sum_ints (nrej_omb_std(:))
209 call da_proc_sum_ints (nrej(:))
210 call da_proc_sum_ints (ngood(:))
213 if (num_fgat_time > 1) then
214 write(filename,'(i2.2,a,i2.2)') it,'_qcstat_'//trim(iv%instid(i)%rttovid_string)//'_',iv%time
216 write(filename,'(i2.2,a)') it,'_qcstat_'//trim(iv%instid(i)%rttovid_string)
219 call da_get_unit(fgat_rad_unit)
220 open(fgat_rad_unit,file=trim(filename),form='formatted',iostat=ios)
222 write(unit=message(1),fmt='(A,A)') 'error opening the output file ', filename
223 call da_error(__FILE__,__LINE__,message(1:1))
226 write(fgat_rad_unit, fmt='(/a/)') ' Quality Control Statistics for '//iv%instid(i)%rttovid_string
227 write(fgat_rad_unit,'(a20,i7)') ' num_proc_domain = ', num_proc_domain
228 write(fgat_rad_unit,'(a20,i7)') ' nrej_mixsurface = ', nrej_mixsurface
229 write(fgat_rad_unit,'(a20,i7)') ' nrej_windowchanl = ', nrej_windowchanl
230 write(fgat_rad_unit,'(a20,i7)') ' nrej_si = ', nrej_si
231 write(fgat_rad_unit,'(a20,i7)') ' nrej_clw = ', nrej_clw
232 write(fgat_rad_unit,'(a20,i7)') ' nrej_topo = ', nrej_topo
233 write(fgat_rad_unit,'(a20,i7)') ' nrej_limb = ', nrej_limb
234 write(fgat_rad_unit,'(a20)') ' nrej_omb_abs(:) = '
235 write(fgat_rad_unit,'(10i7)') nrej_omb_abs(:)
236 write(fgat_rad_unit,'(a20)') ' nrej_omb_std(:) = '
237 write(fgat_rad_unit,'(10i7)') nrej_omb_std(:)
238 write(fgat_rad_unit,'(a20)') ' nrej(:) = '
239 write(fgat_rad_unit,'(10i7)') nrej(:)
240 write(fgat_rad_unit,'(a20)') ' ngood(:) = '
241 write(fgat_rad_unit,'(10i7)') ngood(:)
244 call da_free_unit(fgat_rad_unit)
246 if (trace_use) call da_trace_exit("da_qc_mwts")
248 end subroutine da_qc_mwts