Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_qscat / da_check_max_iv_qscat.inc
blob44b1337fb249533c88662eddee943cd8f8a32454
1 subroutine da_check_max_iv_qscat(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
14    integer,       intent(inout) :: num_qcstat_conv(:,:,:,:)
16    logical :: failed,failed1,failed2
17    integer :: n
18    
19    if (trace_use_dull) call da_trace_entry("da_check_max_iv_qscat")
21    !---------------------------------------------------------------------------
22    ! [1.0] Perform maximum innovation vector check:
23    !---------------------------------------------------------------------------
25    do n=iv%info(qscat)%n1,iv%info(qscat)%n2
26       if(.not. qc_rej_both)then
27          if(wind_sd_qscat)then
28             failed=.false.
29             if( iv%qscat(n)%u%qc >= obs_qc_pointer ) then
30                 call da_max_error_qc (it,iv%info(qscat), n, iv%qscat(n)%u, max_error_spd,failed)
31                 if( iv%info(qscat)%proc_domain(1,n) ) then
32                     num_qcstat_conv(1,qscat,1,1) = num_qcstat_conv(1,qscat,1,1) + 1
33                     if(failed) then
34                        num_qcstat_conv(2,qscat,1,1) = num_qcstat_conv(2,qscat,1,1) + 1
35                        if ( write_rej_obs_conv ) then
36                        write(qcstat_conv_unit,'(2x,a10,2x,a4,2f12.2,a12)')&
37                        'qscat',ob_vars(1),iv%info(qscat)%lat(1,n),iv%info(qscat)%lon(1,n),'1013.25'
38                        end if
39                     end if
40                 end if
41             end if
43             failed=.false.
44             if( iv%qscat(n)%v%qc >= obs_qc_pointer ) then
45                 call da_max_error_qc (it,iv%info(qscat), n, iv%qscat(n)%v, max_error_dir,failed)
46                 if( iv%info(qscat)%proc_domain(1,n) ) then
47                     num_qcstat_conv(1,qscat,2,1) = num_qcstat_conv(1,qscat,2,1) + 1
48                     if(failed)then
49                        num_qcstat_conv(2,qscat,2,1) = num_qcstat_conv(2,qscat,2,1) + 1
50                        if ( write_rej_obs_conv ) then
51                        write(qcstat_conv_unit,'(2x,a10,2x,a4,2f12.2,a12)')&
52                              'qscat',ob_vars(2),iv%info(qscat)%lat(1,n),iv%info(qscat)%lon(1,n),'1013.25'
53                        end if
54                     end if
55                 end if
56             end if
58          else
59             failed=.false.
60             if( iv%qscat(n)%u%qc >= obs_qc_pointer ) then
61                 call da_max_error_qc (it,iv%info(qscat), n, iv%qscat(n)%u, max_error_uv,failed)
62                 if( iv%info(qscat)%proc_domain(1,n) ) then
63                     num_qcstat_conv(1,qscat,1,1) = num_qcstat_conv(1,qscat,1,1) + 1
64                     if(failed) then
65                        num_qcstat_conv(2,qscat,1,1) = num_qcstat_conv(2,qscat,1,1) + 1
66                        if ( write_rej_obs_conv ) then
67                        write(qcstat_conv_unit,'(2x,a10,2x,a4,2f12.2,a12)')&
68                              'qscat',ob_vars(1),iv%info(qscat)%lat(1,n),iv%info(qscat)%lon(1,n),'1013.25'
69                        end if
70                     end if
71                 end if
72             end if
73             failed=.false.
74             if( iv%qscat(n)%v%qc >= obs_qc_pointer ) then
75                 call da_max_error_qc (it,iv%info(qscat), n, iv%qscat(n)%v, max_error_uv,failed)
76                 if( iv%info(qscat)%proc_domain(1,n) ) then
77                     num_qcstat_conv(1,qscat,2,1) = num_qcstat_conv(1,qscat,2,1) + 1
78                     if(failed)then
79                        num_qcstat_conv(2,qscat,2,1) = num_qcstat_conv(2,qscat,2,1) + 1
80                        if ( write_rej_obs_conv ) then
81                        write(qcstat_conv_unit,'(2x,a10,2x,a4,2f12.2,a12)')&
82                              'qscat',ob_vars(2),iv%info(qscat)%lat(1,n),iv%info(qscat)%lon(1,n),'1013.25'
83                        end if
84                     end if
85                 end if
86              end if
87           end if
89           if(wind_sd_qscat)then
90              if(iv%qscat(n)%u%qc == fails_error_max .or. abs(iv%qscat(n)%u%inv) >= max_omb_spd) then
91                 iv%qscat(n)%u%qc = fails_error_max
92                 iv%qscat(n)%u%inv = 0.0
93              endif
94              if(iv%qscat(n)%v%qc == fails_error_max .or. abs(iv%qscat(n)%v%inv) >= max_omb_dir) then
95                 iv%qscat(n)%v%qc = fails_error_max
96                 iv%qscat(n)%v%inv = 0.0
97              endif
98           endif
100        else
102           failed1=.false.
103           failed2=.false.
105           if( iv%qscat(n)%v%qc >= obs_qc_pointer .or. iv%qscat(n)%u%qc >= obs_qc_pointer )  then
106               if(wind_sd_qscat)then
107                  call da_max_error_qc (it,iv%info(qscat), n, iv%qscat(n)%u, max_error_spd,failed1)
108                  call da_max_error_qc (it,iv%info(qscat), n, iv%qscat(n)%v, max_error_dir,failed2)
109               else
110                  call da_max_error_qc (it,iv%info(qscat), n, iv%qscat(n)%u, max_error_uv,failed1)
111                  call da_max_error_qc (it,iv%info(qscat), n, iv%qscat(n)%v, max_error_uv,failed2)
112               endif
114               if( iv%info(qscat)%proc_domain(1,n) ) then
115                   num_qcstat_conv(1,qscat,1,1) = num_qcstat_conv(1,qscat,1,1) + 1
116                   num_qcstat_conv(1,qscat,2,1) = num_qcstat_conv(1,qscat,2,1) + 1
118                   if(failed1 .or. failed2) then
119                      num_qcstat_conv(2,qscat,1,1) = num_qcstat_conv(2,qscat,1,1) + 1
120                      if ( write_rej_obs_conv ) then
121                      write(qcstat_conv_unit,'(2x,a10,2x,a4,2f12.2,a12)')&
122                            'qscat',ob_vars(1),iv%info(qscat)%lat(1,n),iv%info(qscat)%lon(1,n),'1013.25'
123                      end if
124                      num_qcstat_conv(2,qscat,2,1) = num_qcstat_conv(2,qscat,2,1) + 1
125                      if ( write_rej_obs_conv ) then
126                      write(qcstat_conv_unit,'(2x,a10,2x,a4,2f12.2,a12)')&
127                            'qscat',ob_vars(2),iv%info(qscat)%lat(1,n),iv%info(qscat)%lon(1,n),'1013.25'
128                      end if
129                   end if
130               end if
131           end if
133           if(wind_sd_qscat)then
134              if(iv%qscat(n)%u%qc == fails_error_max .or. iv%qscat(n)%v%qc == fails_error_max .or. &
135                 abs(iv%qscat(n)%v%inv) >= max_omb_dir .or. abs(iv%qscat(n)%u%inv) >= max_omb_spd )then
136                 iv%qscat(n)%u%qc = fails_error_max
137                 iv%qscat(n)%v%qc = fails_error_max
138                 iv%qscat(n)%u%inv = 0.0
139                 iv%qscat(n)%v%inv = 0.0
140              endif
141           else
142              if(iv%qscat(n)%u%qc == fails_error_max .or. iv%qscat(n)%v%qc == fails_error_max ) then
143                 iv%qscat(n)%u%qc = fails_error_max
144                 iv%qscat(n)%v%qc = fails_error_max
145                 iv%qscat(n)%u%inv = 0.0
146                 iv%qscat(n)%v%inv = 0.0
147              endif
148          endif
149       endif
151     end do
153    if (trace_use_dull) call da_trace_exit("da_check_max_iv_qscat")
155 end subroutine da_check_max_iv_qscat