1 subroutine da_check_max_iv_tamdar(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(:,:,:,:)
19 if (trace_use_dull) call da_trace_entry("da_check_max_iv_tamdar")
21 !---------------------------------------------------------------------------
22 ! [1.0] Perform maximum innovation vector check:
23 !---------------------------------------------------------------------------
25 do n = iv%info(tamdar)%n1,iv%info(tamdar)%n2
26 do k = 1, iv%info(tamdar)%levels(n)
27 call da_get_print_lvl(iv%tamdar(n)%p(k),ipr)
30 if( iv%tamdar(n)%u(k)%qc >= obs_qc_pointer ) &
31 call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%u(k), max_error_uv,failed)
32 if( iv%info(tamdar)%proc_domain(k,n) ) then
33 num_qcstat_conv(1,tamdar,1,ipr) = num_qcstat_conv(1,tamdar,1,ipr) + 1
35 num_qcstat_conv(2,tamdar,1,ipr) = num_qcstat_conv(2,tamdar,1,ipr) + 1
36 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
37 'tamdar',ob_vars(1),iv%info(tamdar)%lat(k,n),iv%info(tamdar)%lon(k,n),0.01*iv%tamdar(n)%p(k)
42 if( iv%tamdar(n)%v(k)%qc >= obs_qc_pointer ) &
43 call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%v(k), max_error_uv,failed)
44 if( iv%info(tamdar)%proc_domain(k,n) ) then
45 num_qcstat_conv(1,tamdar,2,ipr) = num_qcstat_conv(1,tamdar,2,ipr) + 1
47 num_qcstat_conv(2,tamdar,2,ipr) = num_qcstat_conv(2,tamdar,2,ipr) + 1
48 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
49 'tamdar',ob_vars(2),iv%info(tamdar)%lat(k,n),iv%info(tamdar)%lon(k,n),0.01*iv%tamdar(n)%p(k)
54 if( iv%tamdar(n)%t(k)%qc >= obs_qc_pointer ) &
55 call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%t(k), max_error_t ,failed)
56 if( iv%info(tamdar)%proc_domain(k,n) ) then
57 num_qcstat_conv(1,tamdar,3,ipr) = num_qcstat_conv(1,tamdar,3,ipr) + 1
59 num_qcstat_conv(2,tamdar,3,ipr) = num_qcstat_conv(2,tamdar,3,ipr) + 1
60 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
61 'tamdar',ob_vars(3),iv%info(tamdar)%lat(k,n),iv%info(tamdar)%lon(k,n),0.01*iv%tamdar(n)%p(k)
66 if( iv%tamdar(n)%q(k)%qc >= obs_qc_pointer ) &
67 call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%q(k), max_error_q ,failed)
68 if( iv%info(tamdar)%proc_domain(k,n) ) then
69 num_qcstat_conv(1,tamdar,4,ipr) = num_qcstat_conv(1,tamdar,4,ipr) + 1
71 num_qcstat_conv(2,tamdar,4,ipr) = num_qcstat_conv(2,tamdar,4,ipr) + 1
72 write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
73 'tamdar',ob_vars(4),iv%info(tamdar)%lat(k,n),iv%info(tamdar)%lon(k,n),0.01*iv%tamdar(n)%p(k)
80 if (trace_use_dull) call da_trace_exit("da_check_max_iv_tamdar")
82 end subroutine da_check_max_iv_tamdar