Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_pilot / da_check_max_iv_pilot.inc
blob3a63ab3b5cbc929f4e99cc79aed97762f734e191
1 subroutine da_check_max_iv_pilot(iv, it, num_qcstat_conv)     
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    !-----------------------------------------------------------------------
10    implicit none
12    type(iv_type), intent(inout) :: iv
13    integer,       intent(in)    :: it      ! Outer iteration
14    integer,       intent(inout) :: num_qcstat_conv(:,:,:,:)
16    integer :: k,n,ipr
17    logical :: failed,failed1,failed2
18    
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
30             if(wind_sd_pilot)then
31                failed=.false.
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
36                        if(failed) then
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)
41                           end if
42                        end if
43                    end if
44                 end if
46                 failed=.false.
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
51                         if(failed)then
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)
56                            end if
57                         end if
58                     end if
59                 end if
60              else
61                 failed=.false.
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
66                         if(failed) then
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)
71                            end if
72                         end if
73                     end if
74                 end if
76                 failed=.false.
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
81                         if(failed)then
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)
86                            end if
87                         end if
88                     end if
89                 end if
90               end if
92               if(wind_sd_pilot)then
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
96               endif
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
100               endif
101            endif
103         else
104            failed1=.false.
105            failed2=.false.
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)
111                else
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)
114                endif
115            endif
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)
126                   end if
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)
131                   end if
132                endif
133            endif
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
142               endif
143            else
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
149               endif
150            endif
151        endif
153       end do
154    end do
155    
156    if (trace_use_dull) call da_trace_exit("da_check_max_iv_pilot")
158 end subroutine da_check_max_iv_pilot