da_qc_amsua.inc

References to this file elsewhere.
1 subroutine da_qc_amsua (i, nchan, ob, iv)
2 
3    !---------------------------------------------------------------------------
4    ! Purpose: perform quality control for amsua 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 
15    ! local variables
16    integer   :: n,scanpos,k,isflg,ios,fgat_rad_unit
17    logical   :: lmix
18    real      :: si
19    ! real    :: satzen
20    integer   :: ngood(nchan),nrej(nchan),nrej_omb_abs(nchan), &
21                 nrej_omb_std(nchan),      &
22                 nrej_mixsurface,nrej_windowchanl, nrej_si,    &
23                 nrej_clw,nrej_topo, num_proc_domain,  &
24                 nrej_limb
25 
26    character(len=30)  :: filename
27 
28    if (trace_use) call da_trace_entry("da_qc_amsua")
29 
30    ngood(:)        = 0
31    nrej(:)         = 0
32    nrej_omb_abs(:) = 0
33    nrej_omb_std(:) = 0
34    nrej_mixsurface = 0
35    nrej_windowchanl= 0
36    nrej_si         = 0
37    nrej_clw        = 0
38    nrej_topo       = 0
39    nrej_limb       = 0
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 by flags 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          !  b.  reject channels 1~4 over land/sea-ice/snow
61          !------------------------------------------------------
62          if (isflg > 0) then 
63             iv%instid(i)%tb_qc(1:4,n)  = qc_bad
64             if (iv%instid(i)%info%proc_domain(1,n)) &
65                nrej_windowchanl = nrej_windowchanl + 1
66            ! reject whole pixel if not over sea for global case
67             if (global) iv%instid(i)%tb_qc(:,n)  = qc_bad
68          end if
69 
70          !  c.  reject channels 13,14(above top model 10mb),15 
71          !------------------------------------------------------
72          iv%instid(i)%tb_qc(13:15,n)  = qc_bad
73 
74          !    reject limb obs 
75          !------------------------------------------------------
76          scanpos = iv%instid(i)%scanpos(n)
77          if (scanpos <= 3 .or. scanpos >= 28) then
78             iv%instid(i)%tb_qc(:,n)  =  qc_bad
79             if (iv%instid(i)%info%proc_domain(1,n)) &
80                   nrej_limb = nrej_limb + 1
81          end if
82 
83          ! satzen  = rad%satzen
84          ! if (abs(satzen) > 45.0) iv%instid(i)%tb_qc(:,n)  =  qc_bad
85 
86          !  d. check precipitation 
87          !-----------------------------------------------------------
88          if (ob%instid(i)%tb(1,n) > 0.0 .and. &
89              ob%instid(i)%tb(15,n) > 0.0) then
90             si = ob%instid(i)%tb(1,n) - ob%instid(i)%tb(15,n)
91             if (si >= 3.0) then
92                iv%instid(i)%tb_qc(:,n) = qc_bad
93                iv%instid(i)%cloud_flag(:,n) = qc_bad
94                if (iv%instid(i)%info%proc_domain(1,n)) &
95                   nrej_si = nrej_si + 1
96             end if
97          end if
98 
99          if (iv%instid(i)%clwp(n) >= 0.2) then
100             iv%instid(i)%tb_qc(:,n) = qc_bad
101             iv%instid(i)%cloud_flag(:,n) = qc_bad
102             if (iv%instid(i)%info%proc_domain(1,n)) &
103                nrej_clw = nrej_clw + 1
104          end if
105 
106          !   3.1 Estimate Cloud Liquid Water (CLW) in mm over sea
107          !       (Grody etal. 2001, JGR, Equation 5b,7c,7d,9)
108          !---------------------------------------------------------
109          ! if (isflg == 0) then
110          !    coszen =  cos(iv%instid(i)%satzen(n))
111          !    d0     =  8.24-(2.622-1.846*coszen)*coszen
112          !    d1     =  0.754
113          !    d2     =  -2.265
114          !    ts     =  iv%instid(i)%ts(n)
115          !    tb1    =  ob%instid(i)%tb(1,n)
116          !    tb2    =  ob%instid(i)%tb(2,n)
117          !    clw    =  coszen*(d0+d1*log(ts-tb1)+d2*log(ts-tb2))
118          !    clw    =  clw - 0.03
119          ! end if
120 
121 
122          !  e. check surface height/pressure
123          !-----------------------------------------------------------
124          ! sfchgt = ivrad%info(n)%elv
125          ! if (sfchgt >=) then
126          ! else 
127          ! end if
128 
129          if ((isflg .ne. 0) .and. (iv%instid(i)%ps(n) < 850.0)) then
130             iv%instid(i)%tb_qc(5,n)  = qc_bad
131             if (iv%instid(i)%info%proc_domain(1,n)) &
132                nrej_topo = nrej_topo + 1
133          end if
134 
135          !  g. check iuse
136          !-----------------------------------------------------------
137          do k = 1, nchan
138             if (satinfo(i)%iuse(k) .eq. -1) &
139                iv%instid(i)%tb_qc(k,n)  = qc_bad
140          end do
141 
142          !  f. check innovation
143          !-----------------------------------------------------------
144          do k = 1, nchan
145 
146          ! absolute departure check
147             if (abs(iv%instid(i)%tb_inv(k,n)) > 15.0) then
148                iv%instid(i)%tb_qc(k,n)  = qc_bad
149                if (iv%instid(i)%info%proc_domain(1,n)) &
150                   nrej_omb_abs(k) = nrej_omb_abs(k) + 1
151             end if
152 
153          ! relative departure check
154             if (use_error_factor_rad) then
155                iv%instid(i)%tb_error(k,n) = &
156                    satinfo(i)%error_std(k)*satinfo(i)%error_factor(k)
157             else
158                iv%instid(i)%tb_error(k,n) = satinfo(i)%error_std(k)
159             end if
160 
161             if (abs(iv%instid(i)%tb_inv(k,n)) > 3.0*iv%instid(i)%tb_error(k,n)) then
162                 iv%instid(i)%tb_qc(k,n)  = qc_bad
163                 if (iv%instid(i)%info%proc_domain(1,n)) &
164                    nrej_omb_std(k) = nrej_omb_std(k) + 1
165             end if
166 
167          ! final QC decsion
168             if (iv%instid(i)%tb_qc(k,n) == qc_bad) then
169                iv%instid(i)%tb_error(k,n) = 500.0
170                if (iv%instid(i)%info%proc_domain(1,n)) &
171                   nrej(k) = nrej(k) + 1
172             else
173                if (iv%instid(i)%info%proc_domain(1,n)) &
174                   ngood(k) = ngood(k) + 1
175             end if
176          end do ! chan
177       end do ! end loop pixel
178  
179    ! Do inter-processor communication to gather statistics.
180    call da_proc_sum_int (num_proc_domain)
181    call da_proc_sum_int (nrej_mixsurface)
182    call da_proc_sum_int (nrej_windowchanl)
183    call da_proc_sum_int (nrej_si )
184    call da_proc_sum_int (nrej_clw)
185    call da_proc_sum_int (nrej_topo)
186    call da_proc_sum_int (nrej_limb)
187    call da_proc_sum_ints (nrej_omb_abs(:))
188    call da_proc_sum_ints (nrej_omb_std(:))
189    call da_proc_sum_ints (nrej(:))
190    call da_proc_sum_ints (ngood(:))
191 
192    if (rootproc) then
193       if (num_fgat_time > 1) then
194          write(filename,'(a,i2.2)') 'qcstat_'//trim(iv%instid(i)%rttovid_string)//'_',iv%time
195       else
196          filename = 'qcstat_'//trim(iv%instid(i)%rttovid_string)
197       end if
198 
199       call da_get_unit(fgat_rad_unit)
200       open(fgat_rad_unit,file=trim(filename),form='formatted',iostat=ios)
201       if (ios /= 0) then
202          write(unit=message(1),fmt='(A,A)') 'error opening the output file ', filename
203          call da_error(__FILE__,__LINE__,message(1:1))
204       end if
205 
206       write(fgat_rad_unit, fmt='(/a/)') ' Quality Control Statistics for '//iv%instid(i)%rttovid_string
207       write(fgat_rad_unit,'(a20,i7)') ' num_proc_domain  = ', num_proc_domain
208       write(fgat_rad_unit,'(a20,i7)') ' nrej_mixsurface  = ', nrej_mixsurface
209       write(fgat_rad_unit,'(a20,i7)') ' nrej_windowchanl = ', nrej_windowchanl
210       write(fgat_rad_unit,'(a20,i7)') ' nrej_si          = ', nrej_si
211       write(fgat_rad_unit,'(a20,i7)') ' nrej_clw         = ', nrej_clw
212       write(fgat_rad_unit,'(a20,i7)') ' nrej_topo        = ', nrej_topo
213       write(fgat_rad_unit,'(a20,i7)') ' nrej_limb        = ', nrej_limb
214       write(fgat_rad_unit,'(a20)')    ' nrej_omb_abs(:)  = '
215       write(fgat_rad_unit,'(10i7)')     nrej_omb_abs(:)
216       write(fgat_rad_unit,'(a20)')    ' nrej_omb_std(:)  = '
217       write(fgat_rad_unit,'(10i7)')     nrej_omb_std(:)
218       write(fgat_rad_unit,'(a20)')    ' nrej(:)          = '
219       write(fgat_rad_unit,'(10i7)')     nrej(:)
220       write(fgat_rad_unit,'(a20)')    ' ngood(:)         = '
221       write(fgat_rad_unit,'(10i7)')     ngood(:)
222 
223       close(fgat_rad_unit)
224       call da_free_unit(fgat_rad_unit)
225    end if
226 
227    if (trace_use) call da_trace_exit("da_qc_amsua")
228 
229 end subroutine da_qc_amsua
230 
231