updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_radiance / da_qc_amsub.inc
blob6173607d56523d2ccc4bffa50e5e6d1f6220eb49
1 subroutine da_qc_amsub (it, i, nchan, ob, iv)
3    !---------------------------------------------------------------------------
4    ! Purpose: perform quality control for amsub data.
5    !---------------------------------------------------------------------------
7    implicit none
9    integer, intent(in)             :: it         ! outer loop count
10    integer, intent(in)             :: i          ! sensor index.
11    integer, intent(in)             :: nchan      ! number of channel
12    type (y_type),  intent(in)      :: ob         ! Observation structure.
13    type (iv_type), intent(inout)   :: iv         ! O-B structure.
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
26    character(len=30)  :: filename
28    if (trace_use) call da_trace_entry("da_qc_amsub")
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
41    num_proc_domain = 0
43       do n= iv%instid(i)%info%n1,iv%instid(i)%info%n2
45          if (iv%instid(i)%info%proc_domain(1,n)) &
46                num_proc_domain = num_proc_domain + 1
48          ! 0.0  initialise QC flags by assuming good obs
49          !---------------------------------------------
50          iv%instid(i)%tb_qc(:,n) = qc_good
51          if (crtm_cloud) go to 2508
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 >= 83) 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 .and. &
88               ob%instid(i)%tb(2,n) > 0.0) then
89             si = ob%instid(i)%tb(1,n) - ob%instid(i)%tb(2,n)
90             if (si >= 3.0) then
91                iv%instid(i)%tb_qc(:,n) = qc_bad
92                iv%instid(i)%cloud_flag(:,n) = qc_bad
93                if (iv%instid(i)%info%proc_domain(1,n)) &
94                   nrej_si = nrej_si + 1
95             end if
96          end if
98          if (iv%instid(i)%clwp(n) >= 0.2) then
99             iv%instid(i)%tb_qc(:,n) = qc_bad
100             iv%instid(i)%cloud_flag(:,n) = qc_bad
101             if (iv%instid(i)%info%proc_domain(1,n)) &
102                nrej_clw = nrej_clw + 1
103          end if
105          !  e. check surface height/pressure
106          !------------------------------------------------------
107          ! sfchgt = ivrad%info(n)%elv
108          ! if (sfchgt >=) then
109          !
110          ! else 
111          ! end if
113          if ((isflg .ne. 0) .and. (iv%instid(i)%ps(n) < 800.0)) then
114             iv%instid(i)%tb_qc(5,n)  = qc_bad
115             if (iv%instid(i)%info%proc_domain(1,n)) &
116                nrej_topo = nrej_topo + 1
117          end if
119          !  g. check iuse (pre-rejected channels by .info files)
120          !------------------------------------------------------
121          do k = 1, nchan
122            if (satinfo(i)%iuse(k) .eq. -1) &
123                 iv%instid(i)%tb_qc(k,n) = qc_bad
124          end do
126 2508      continue
128          !  f. check innovation
129          !-----------------------------------------------------------
130          do k = 1, nchan
131           ! absolute departure check
132           if (.not. crtm_cloud) then
133             if (abs(iv%instid(i)%tb_inv(k,n)) > 15.0) then
134                iv%instid(i)%tb_qc(k,n)  = qc_bad
135                if (iv%instid(i)%info%proc_domain(1,n)) &
136                   nrej_omb_abs(k) = nrej_omb_abs(k) + 1
137             end if
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 (.not. crtm_cloud) then
149             if (abs(iv%instid(i)%tb_inv(k,n)) > 3.0*iv%instid(i)%tb_error(k,n)) then
150                iv%instid(i)%tb_qc(k,n)  = qc_bad
151                if (iv%instid(i)%info%proc_domain(1,n)) &
152                      nrej_omb_std(k) = nrej_omb_std(k) + 1
153             end if
155            ! final QC decision
156             if (iv%instid(i)%tb_qc(k,n) == qc_bad) then
157                iv%instid(i)%tb_error(k,n) = 500.0
158                if (iv%instid(i)%info%proc_domain(1,n)) &
159                   nrej(k) = nrej(k) + 1
160             else
161                if (iv%instid(i)%info%proc_domain(1,n)) &
162                   ngood(k) = ngood(k) + 1
163             end if
164           end if
166          end do ! chan
167       end do ! end loop pixel
169    ! Do inter-processor communication to gather statistics.
171    call da_proc_sum_int (num_proc_domain)
172    call da_proc_sum_int (nrej_mixsurface)
173    call da_proc_sum_int (nrej_windowchanl)
174    call da_proc_sum_int (nrej_si)
175    call da_proc_sum_int (nrej_clw)
176    call da_proc_sum_int (nrej_topo)
177    call da_proc_sum_int (nrej_limb)
178    call da_proc_sum_ints (nrej_omb_abs(:))
179    call da_proc_sum_ints (nrej_omb_std(:))
180    call da_proc_sum_ints (nrej(:))
181    call da_proc_sum_ints (ngood(:))
183    if (rootproc) then
184       if (num_fgat_time > 1) then
185          write(filename,'(i2.2,a,i2.2)') it, '_qcstat_'//trim(iv%instid(i)%rttovid_string)//'_',iv%time
186       else
187          write(filename,'(i2.2,a)') it, '_qcstat_'//trim(iv%instid(i)%rttovid_string)
188       end if
189       call da_get_unit(fgat_rad_unit)
190       open(fgat_rad_unit,file=trim(filename),form='formatted',iostat=ios)
191       if (ios /= 0) then
192          write(unit=message(1),fmt='(A,A)') 'error opening the output file ', filename
193          call da_error(__FILE__,__LINE__,message(1:1))
194       end if
196       write(fgat_rad_unit, fmt='(/a/)') &
197          'Quality Control Statistics for '//iv%instid(i)%rttovid_string
198       write(fgat_rad_unit,'(a20,i7)') ' num_proc_domain  = ', num_proc_domain
199       write(fgat_rad_unit,'(a20,i7)') ' nrej_mixsurface  = ', nrej_mixsurface
200       write(fgat_rad_unit,'(a20,i7)') ' nrej_windowchanl = ', nrej_windowchanl
201       write(fgat_rad_unit,'(a20,i7)') ' nrej_si          = ', nrej_si
202       write(fgat_rad_unit,'(a20,i7)') ' nrej_clw         = ', nrej_clw
203       write(fgat_rad_unit,'(a20,i7)') ' nrej_topo        = ', nrej_topo
204       write(fgat_rad_unit,'(a20,i7)') ' nrej_limb        = ', nrej_limb
205       write(fgat_rad_unit,'(a20)')    ' nrej_omb_abs(:)  = '
206       write(fgat_rad_unit,'(10i7)')     nrej_omb_abs(:)
207       write(fgat_rad_unit,'(a20)')    ' nrej_omb_std(:)  = '
208       write(fgat_rad_unit,'(10i7)')     nrej_omb_std(:)
209       write(fgat_rad_unit,'(a20)')    ' nrej(:)          = '
210       write(fgat_rad_unit,'(10i7)')     nrej(:)
211       write(fgat_rad_unit,'(a20)')    ' ngood(:)         = '
212       write(fgat_rad_unit,'(10i7)')     ngood(:)
214       close(fgat_rad_unit)
215       call da_free_unit(fgat_rad_unit)
216    end if
218    if (trace_use) call da_trace_exit("da_qc_amsub")
220 end subroutine da_qc_amsub