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