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