1 subroutine da_check_max_iv_pilot(iv, it, num_qcstat_conv)
3 !-----------------------------------------------------------------------
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 !-----------------------------------------------------------------------
12 type(iv_type), intent(inout) :: iv
13 integer, intent(in) :: it ! Outer iteration
14 integer, intent(inout) :: num_qcstat_conv(:,:,:,:)
17 logical :: failed,failed1,failed2
19 if (trace_use_dull) call da_trace_entry("da_check_max_iv_pilot")
21 !---------------------------------------------------------------------------
22 ! [1.0] Perform maximum innovation vector check:
23 !---------------------------------------------------------------------------
25 do n = iv%info(pilot)%n1,iv%info(pilot)%n2
26 do k = 1, iv%info(pilot)%levels(n)
27 call da_get_print_lvl(iv%pilot(n)%p(k),ipr)
29 if(.not. qc_rej_both)then
32 if( iv%pilot(n)%u(k)%qc >= obs_qc_pointer ) then
33 call da_max_error_qc (it,iv%info(pilot), n, iv%pilot(n)%u(k), max_error_spd,failed)
34 if( iv%info(pilot)%proc_domain(k,n) ) then
35 num_qcstat_conv(1,pilot,1,ipr) = num_qcstat_conv(1,pilot,1,ipr) + 1
37 num_qcstat_conv(2,pilot,1,ipr) = num_qcstat_conv(2,pilot,1,ipr) + 1
38 if ( write_rej_obs_conv ) then
39 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
40 'pilot',ob_vars(1),iv%info(pilot)%lat(k,n),iv%info(pilot)%lon(k,n),0.01*iv%pilot(n)%p(k)
47 if( iv%pilot(n)%v(k)%qc >= obs_qc_pointer ) then
48 call da_max_error_qc (it,iv%info(pilot), n, iv%pilot(n)%v(k), max_error_dir,failed)
49 if( iv%info(pilot)%proc_domain(k,n) ) then
50 num_qcstat_conv(1,pilot,2,ipr) = num_qcstat_conv(1,pilot,2,ipr) + 1
52 num_qcstat_conv(2,pilot,2,ipr) = num_qcstat_conv(2,pilot,2,ipr) + 1
53 if ( write_rej_obs_conv ) then
54 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
55 'pilot',ob_vars(2),iv%info(pilot)%lat(k,n),iv%info(pilot)%lon(k,n),0.01*iv%pilot(n)%p(k)
62 if( iv%pilot(n)%u(k)%qc >= obs_qc_pointer ) then
63 call da_max_error_qc (it,iv%info(pilot), n, iv%pilot(n)%u(k), max_error_uv,failed)
64 if( iv%info(pilot)%proc_domain(k,n) ) then
65 num_qcstat_conv(1,pilot,1,ipr) = num_qcstat_conv(1,pilot,1,ipr) + 1
67 num_qcstat_conv(2,pilot,1,ipr) = num_qcstat_conv(2,pilot,1,ipr) + 1
68 if ( write_rej_obs_conv ) then
69 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
70 'pilot',ob_vars(1),iv%info(pilot)%lat(k,n),iv%info(pilot)%lon(k,n),0.01*iv%pilot(n)%p(k)
77 if( iv%pilot(n)%v(k)%qc >= obs_qc_pointer ) then
78 call da_max_error_qc (it,iv%info(pilot), n, iv%pilot(n)%v(k), max_error_uv,failed)
79 if( iv%info(pilot)%proc_domain(k,n) ) then
80 num_qcstat_conv(1,pilot,2,ipr) = num_qcstat_conv(1,pilot,2,ipr) + 1
82 num_qcstat_conv(2,pilot,2,ipr) = num_qcstat_conv(2,pilot,2,ipr) + 1
83 if ( write_rej_obs_conv ) then
84 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
85 'pilot',ob_vars(2),iv%info(pilot)%lat(k,n),iv%info(pilot)%lon(k,n),0.01*iv%pilot(n)%p(k)
93 if(iv%pilot(n)%u(k)%qc == fails_error_max .or. abs(iv%pilot(n)%u(k)%inv) >= max_omb_spd) then
94 iv%pilot(n)%u(k)%qc = fails_error_max
95 iv%pilot(n)%u(k)%inv = 0.0
97 if(iv%pilot(n)%v(k)%qc == fails_error_max .or. abs(iv%pilot(n)%v(k)%inv) >= max_omb_dir) then
98 iv%pilot(n)%v(k)%qc = fails_error_max
99 iv%pilot(n)%v(k)%inv = 0.0
107 if( iv%pilot(n)%v(k)%qc >= obs_qc_pointer .or. iv%pilot(n)%u(k)%qc >= obs_qc_pointer ) then
108 if(wind_sd_pilot)then
109 call da_max_error_qc (it,iv%info(pilot), n, iv%pilot(n)%u(k), max_error_spd,failed1)
110 call da_max_error_qc (it,iv%info(pilot), n, iv%pilot(n)%v(k), max_error_dir,failed2)
112 call da_max_error_qc (it,iv%info(pilot), n, iv%pilot(n)%u(k), max_error_uv,failed1)
113 call da_max_error_qc (it,iv%info(pilot), n, iv%pilot(n)%v(k), max_error_uv,failed2)
117 if( iv%info(pilot)%proc_domain(k,n) ) then
118 num_qcstat_conv(1,pilot,1,ipr) = num_qcstat_conv(1,pilot,1,ipr) + 1
119 num_qcstat_conv(1,pilot,2,ipr) = num_qcstat_conv(1,pilot,2,ipr) + 1
121 if(failed1 .or. failed2) then
122 num_qcstat_conv(2,pilot,1,ipr) = num_qcstat_conv(2,pilot,1,ipr) + 1
123 if ( write_rej_obs_conv ) then
124 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
125 'pilot',ob_vars(1),iv%info(pilot)%lat(k,n),iv%info(pilot)%lon(k,n),0.01*iv%pilot(n)%p(k)
127 num_qcstat_conv(2,pilot,2,ipr) = num_qcstat_conv(2,pilot,2,ipr) + 1
128 if ( write_rej_obs_conv ) then
129 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
130 'pilot',ob_vars(2),iv%info(pilot)%lat(k,n),iv%info(pilot)%lon(k,n),0.01*iv%pilot(n)%p(k)
135 if(wind_sd_pilot)then
136 if(iv%pilot(n)%u(k)%qc == fails_error_max .or. iv%pilot(n)%v(k)%qc == fails_error_max .or. &
137 abs(iv%pilot(n)%v(k)%inv) >= max_omb_dir .or. abs(iv%pilot(n)%u(k)%inv) >= max_omb_spd )then
138 iv%pilot(n)%u(k)%qc = fails_error_max
139 iv%pilot(n)%v(k)%qc = fails_error_max
140 iv%pilot(n)%u(k)%inv = 0.0
141 iv%pilot(n)%v(k)%inv = 0.0
144 if(iv%pilot(n)%u(k)%qc == fails_error_max .or. iv%pilot(n)%v(k)%qc == fails_error_max ) then
145 iv%pilot(n)%u(k)%qc = fails_error_max
146 iv%pilot(n)%v(k)%qc = fails_error_max
147 iv%pilot(n)%u(k)%inv = 0.0
148 iv%pilot(n)%v(k)%inv = 0.0
156 if (trace_use_dull) call da_trace_exit("da_check_max_iv_pilot")
158 end subroutine da_check_max_iv_pilot