updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_gpsref / da_check_max_iv_gpsref.inc
blob657a9fa30957adc9b4f5d41146c691298b9f7541
1 subroutine da_check_max_iv_gpsref(iv,it, num_qcstat_conv, opt)        
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    ! Update:
6    !    Removed Outerloop check as it is done in da_get_innov
7    !    Author: Syed RH Rizvi,  MMM/NESL/NCAR,  Date: 07/12/2009
8    !
9    !    Added argument: opt = 1, max error check, and count the total num. of obs;
10    !                    opt = 2, only count the rejected num. of obs by all of QC
11    !                             procedures specified for GPSREF.
12    !        
13    !                             Shu-Ya Chen and Y.-R. Guo, 12/30/2011  
14    !-----------------------------------------------------------------------
16    implicit none
18    type(iv_type), intent(inout) :: iv
19    integer,       intent(in)    :: it       ! External iteration.
20    integer,       intent(inout) :: num_qcstat_conv(:,:,:,:)  
22    integer                           :: k, n, ipr
23    logical                           :: failed
24    integer, intent(in) :: opt ! 1: only counting, 2: print out  (sychen add)
25    integer, parameter :: qc_below = -31, qc_middle = -32, qc_above = -33
26    integer, parameter :: qc_step1 = -34, qc_step2  = -35  ! refer to Poli et al. (2009)
27    integer, parameter :: qc_cutoff = -36
29    if (trace_use_dull) call da_trace_entry("da_check_max_iv_gpsref")
31    !----------------------------------------------------------------------------
32    ! [1.0] Perform maximum innovation vector check:
33    !----------------------------------------------------------------------------
34          IF ( opt == 1 ) THEN
35    do n = iv%info(gpsref)%n1,iv%info(gpsref)%n2
36       do k = 1, iv%info(gpsref)%levels(n)
37         call da_get_print_lvl(iv%gpsref(n)%p(k)%inv,ipr)
38         if (iv%gpsref(n)%p(k)%inv == missing_r) ipr = 1
39         failed=.false.
40         if( iv%gpsref(n)%ref(k)%qc >= obs_qc_pointer ) &
41         call da_max_error_qc(it, iv%info(gpsref), n, iv%gpsref(n)%ref(k), max_error_ref, failed)  
42         if( iv%info(gpsref)%proc_domain(k,n) ) &
43                  num_qcstat_conv(1,gpsref,8,ipr) = num_qcstat_conv(1,gpsref,8,ipr) + 1
44       end do
45    end do
46          ENDIF
47   
48          IF ( opt == 2 ) THEN
49    do n = iv%info(gpsref)%n1,iv%info(gpsref)%n2
50       do k = 1, iv%info(gpsref)%levels(n)
51         call da_get_print_lvl(iv%gpsref(n)%p(k)%inv,ipr)
52         if (iv%gpsref(n)%p(k)%inv == missing_r) ipr = 1
53         failed=.false.
54         if ( ( iv%gpsref(n)%ref(k)%qc == fails_error_max ) .or. &
55              ( iv%gpsref(n)%ref(k)%qc == qc_below ) .or. &
56              ( iv%gpsref(n)%ref(k)%qc == qc_middle ) .or. &
57              ( iv%gpsref(n)%ref(k)%qc == qc_above ) .or. &
58              ( iv%gpsref(n)%ref(k)%qc == qc_step1 ) .or. &
59              ( iv%gpsref(n)%ref(k)%qc == qc_step2 ) .or. &
60              ( iv%gpsref(n)%ref(k)%qc == qc_cutoff ) .or. & !hcl-202006
61              ( iv%gpsref(n)%ref(k)%qc == missing_data ) ) then 
62              failed=.true.
63         endif
64             if(failed) then
65             num_qcstat_conv(2,gpsref,8,ipr) = num_qcstat_conv(2,gpsref,8,ipr) + 1
66             if ( write_rej_obs_conv ) then
67             write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2,I5)')&
68              'gpsref',ob_vars(8),iv%info(gpsref)%lat(k,n), &
69              iv%info(gpsref)%lon(k,n),0.01*iv%gpsref(n)%p(k)%inv, &
70              iv%gpsref(n)%ref(k)%qc 
71             end if
72             end if
73       end do
74    end do
75          ENDIF
77    if (trace_use_dull) call da_trace_exit("da_check_max_iv_gpsref")
79 end subroutine da_check_max_iv_gpsref