1 subroutine da_check_max_iv_bogus(iv,ob, it, num_qcstat_conv)
3 !-----------------------------------------------------------------------
5 ! Removed Outerloop check as it is done in da_get_innov
6 ! Author: Syed RH Rizvi, MMM/NESL/NCAR, Date: 07/12/2009
7 !-----------------------------------------------------------------------
11 type(iv_type), intent(inout) :: iv
12 integer, intent(in) :: it ! Outer iteration
13 integer, intent(inout) :: num_qcstat_conv(:,:,:,:)
14 type(y_type), intent(in) :: ob ! Observation structure.
21 if (trace_use_dull) call da_trace_entry("da_check_max_iv_bogus")
23 !---------------------------------------------------------------------------
24 ! [1.0] Perform maximum innovation vector check:
25 !---------------------------------------------------------------------------
27 do n = iv%info(bogus)%n1,iv%info(bogus)%n2
28 do k = 1, iv%info(bogus)%levels(n)
29 call da_get_print_lvl(iv%bogus(n)%p(k),ipr)
32 if( iv%bogus(n)%u(k)%qc >= obs_qc_pointer ) then
33 call da_max_error_qc (it,iv%info(bogus), n, iv%bogus(n)%u(k), max_error_buv,failed)
34 if( iv%info(bogus)%proc_domain(k,n) ) then
35 num_qcstat_conv(1,bogus,1,ipr) = num_qcstat_conv(1,bogus,1,ipr) + 1
37 num_qcstat_conv(2,bogus,1,ipr) = num_qcstat_conv(2,bogus,1,ipr) + 1
38 if ( write_rej_obs_conv ) then
39 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
40 'bogus',ob_vars(1),iv%info(bogus)%lat(k,n),iv%info(bogus)%lon(k,n),0.01*iv%bogus(n)%p(k)
47 if( iv%bogus(n)%v(k)%qc >= obs_qc_pointer ) then
48 call da_max_error_qc (it,iv%info(bogus), n, iv%bogus(n)%v(k), max_error_buv,failed)
49 if( iv%info(bogus)%proc_domain(k,n) ) then
50 num_qcstat_conv(1,bogus,2,ipr) = num_qcstat_conv(1,bogus,2,ipr) + 1
52 num_qcstat_conv(2,bogus,2,ipr) = num_qcstat_conv(2,bogus,2,ipr) + 1
53 if ( write_rej_obs_conv ) then
54 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
55 'bogus',ob_vars(2),iv%info(bogus)%lat(k,n),iv%info(bogus)%lon(k,n),0.01*iv%bogus(n)%p(k)
62 if( iv%bogus(n)%t(k)%qc >= obs_qc_pointer ) then
63 call da_max_error_qc (it,iv%info(bogus), n, iv%bogus(n)%t(k), max_error_bt ,failed)
64 if( iv%info(bogus)%proc_domain(k,n) ) then
65 num_qcstat_conv(1,bogus,3,ipr) = num_qcstat_conv(1,bogus,3,ipr) + 1
67 num_qcstat_conv(2,bogus,3,ipr) = num_qcstat_conv(2,bogus,3,ipr) + 1
68 if ( write_rej_obs_conv ) then
69 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
70 'bogus',ob_vars(3),iv%info(bogus)%lat(k,n),iv%info(bogus)%lon(k,n),0.01*iv%bogus(n)%p(k)
77 if( iv%bogus(n)%q(k)%qc >= obs_qc_pointer ) then
78 call da_max_error_qc (it,iv%info(bogus), n, iv%bogus(n)%q(k), max_error_bq ,failed)
79 if( iv%info(bogus)%proc_domain(k,n) ) then
80 num_qcstat_conv(1,bogus,4,ipr) = num_qcstat_conv(1,bogus,4,ipr) + 1
82 num_qcstat_conv(2,bogus,4,ipr) = num_qcstat_conv(2,bogus,4,ipr) + 1
83 if ( write_rej_obs_conv ) then
84 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
85 'bogus',ob_vars(4),iv%info(bogus)%lat(k,n),iv%info(bogus)%lon(k,n),0.01*iv%bogus(n)%p(k)
94 if( iv%info(bogus)%proc_domain(1,n) ) then
96 if( iv%bogus(n)%slp%qc >= obs_qc_pointer ) then
97 call da_max_error_qc (it,iv%info(bogus), n, iv%bogus(n)%slp, max_error_slp ,failed)
98 num_qcstat_conv(1,bogus,5,1) = num_qcstat_conv(1,bogus,5,1) + 1
100 num_qcstat_conv(2,bogus,5,1) = num_qcstat_conv(2,bogus,5,1) + 1
101 if ( write_rej_obs_conv ) then
102 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
103 'bogus',ob_vars(5),iv%info(bogus)%lat(1,n),iv%info(bogus)%lon(1,n),ob%bogus(n)%slp
111 if (trace_use_dull) call da_trace_exit("da_check_max_iv_bogus")
113 end subroutine da_check_max_iv_bogus