Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / var / da / da_tamdar / da_check_max_iv_tamdar.inc
blob301b69e62e25fb3eb78149d143203acee30d7841
1 subroutine da_check_max_iv_tamdar(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
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) 
29          if(.not. qc_rej_both)then
30             if(wind_sd_tamdar)then
31                failed=.false.
32                if( iv%tamdar(n)%u(k)%qc >= obs_qc_pointer ) then
33                    call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%u(k), max_error_spd,failed)
34                    if( iv%info(tamdar)%proc_domain(k,n) ) then
35                        num_qcstat_conv(1,tamdar,1,ipr) = num_qcstat_conv(1,tamdar,1,ipr) + 1
36                        if(failed) then
37                           num_qcstat_conv(2,tamdar,1,ipr) = num_qcstat_conv(2,tamdar,1,ipr) + 1
38                           if ( write_rej_obs_conv ) then
39                           write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
40                           'tamdar',ob_vars(1),iv%info(tamdar)%lat(k,n),iv%info(tamdar)%lon(k,n),0.01*iv%tamdar(n)%p(k)
41                           end if
42                        end if
43                    end if
44                 end if
46                 failed=.false.
47                 if( iv%tamdar(n)%v(k)%qc >= obs_qc_pointer ) then
48                     call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%v(k), max_error_dir,failed)
49                     if( iv%info(tamdar)%proc_domain(k,n) ) then
50                         num_qcstat_conv(1,tamdar,2,ipr) = num_qcstat_conv(1,tamdar,2,ipr) + 1
51                         if(failed)then
52                            num_qcstat_conv(2,tamdar,2,ipr) = num_qcstat_conv(2,tamdar,2,ipr) + 1
53                            if ( write_rej_obs_conv ) then
54                            write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
55                            'tamdar',ob_vars(2),iv%info(tamdar)%lat(k,n),iv%info(tamdar)%lon(k,n),0.01*iv%tamdar(n)%p(k)
56                            end if
57                         end if
58                     end if
59                 end if
60              else
61                 failed=.false.
62                 if( iv%tamdar(n)%u(k)%qc >= obs_qc_pointer ) then
63                     call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%u(k), max_error_uv,failed)
64                     if( iv%info(tamdar)%proc_domain(k,n) ) then
65                         num_qcstat_conv(1,tamdar,1,ipr) = num_qcstat_conv(1,tamdar,1,ipr) + 1
66                         if(failed) then
67                            num_qcstat_conv(2,tamdar,1,ipr) = num_qcstat_conv(2,tamdar,1,ipr) + 1
68                            if ( write_rej_obs_conv ) then
69                            write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
70                            'tamdar',ob_vars(1),iv%info(tamdar)%lat(k,n),iv%info(tamdar)%lon(k,n),0.01*iv%tamdar(n)%p(k)
71                            end if
72                          end if
73                     end if
74                  end if
76                  failed=.false.
77                  if( iv%tamdar(n)%v(k)%qc >= obs_qc_pointer ) then
78                      call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%v(k), max_error_uv,failed)
79                      if( iv%info(tamdar)%proc_domain(k,n) ) then
80                          num_qcstat_conv(1,tamdar,2,ipr) = num_qcstat_conv(1,tamdar,2,ipr) + 1
81                          if(failed)then
82                             num_qcstat_conv(2,tamdar,2,ipr) = num_qcstat_conv(2,tamdar,2,ipr) + 1
83                             if ( write_rej_obs_conv ) then
84                             write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
85                             'tamdar',ob_vars(2),iv%info(tamdar)%lat(k,n),iv%info(tamdar)%lon(k,n),0.01*iv%tamdar(n)%p(k)
86                             end if
87                          end if
88                      end if
89                  end if
90               end if
92               if(wind_sd_tamdar)then
93                  if(iv%tamdar(n)%u(k)%qc == fails_error_max .or. abs(iv%tamdar(n)%u(k)%inv) >= max_omb_spd) then
94                     iv%tamdar(n)%u(k)%qc = fails_error_max
95                     iv%tamdar(n)%u(k)%inv = 0.0
96                  endif
97                  if(iv%tamdar(n)%v(k)%qc == fails_error_max .or. abs(iv%tamdar(n)%v(k)%inv) >= max_omb_dir) then
98                     iv%tamdar(n)%v(k)%qc = fails_error_max
99                     iv%tamdar(n)%v(k)%inv = 0.0
100                  endif
101               endif
102            else
103               failed1=.false.
104               failed2=.false.
106               if( iv%tamdar(n)%v(k)%qc >= obs_qc_pointer .or. iv%tamdar(n)%u(k)%qc >= obs_qc_pointer )  then
107                   if(wind_sd_tamdar)then
108                      call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%u(k), max_error_spd,failed1)
109                      call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%v(k), max_error_dir,failed2)
110                   else
111                      call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%u(k), max_error_uv,failed1)
112                      call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%v(k), max_error_uv,failed2)
113                   endif
114               endif
116               if( iv%info(tamdar)%proc_domain(k,n) ) then
117                   num_qcstat_conv(1,tamdar,1,ipr) = num_qcstat_conv(1,tamdar,1,ipr) + 1
118                   num_qcstat_conv(1,tamdar,2,ipr) = num_qcstat_conv(1,tamdar,2,ipr) + 1
120                   if(failed1 .or. failed2) then
121                      num_qcstat_conv(2,tamdar,1,ipr) = num_qcstat_conv(2,tamdar,1,ipr) + 1
122                      if ( write_rej_obs_conv ) then
123                      write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
124                            'tamdar',ob_vars(1),iv%info(tamdar)%lat(k,n),iv%info(tamdar)%lon(k,n),0.01*iv%tamdar(n)%p(k)
125                      end if
126                      num_qcstat_conv(2,tamdar,2,ipr) = num_qcstat_conv(2,tamdar,2,ipr) + 1
127                      if ( write_rej_obs_conv ) then
128                      write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
129                            'tamdar',ob_vars(2),iv%info(tamdar)%lat(k,n),iv%info(tamdar)%lon(k,n),0.01*iv%tamdar(n)%p(k)
130                      end if
131                   endif
132                endif
134                if(wind_sd_tamdar)then
135                   if(iv%tamdar(n)%u(k)%qc == fails_error_max .or. iv%tamdar(n)%v(k)%qc == fails_error_max .or. &
136                      abs(iv%tamdar(n)%v(k)%inv) >= max_omb_dir .or. abs(iv%tamdar(n)%u(k)%inv) >= max_omb_spd )then
137                      iv%tamdar(n)%u(k)%qc = fails_error_max
138                      iv%tamdar(n)%v(k)%qc = fails_error_max
139                      iv%tamdar(n)%u(k)%inv = 0.0
140                      iv%tamdar(n)%v(k)%inv = 0.0
141                   endif
142                else
143                   if(iv%tamdar(n)%u(k)%qc == fails_error_max .or. iv%tamdar(n)%v(k)%qc == fails_error_max ) then
144                      iv%tamdar(n)%u(k)%qc = fails_error_max
145                      iv%tamdar(n)%v(k)%qc = fails_error_max
146                      iv%tamdar(n)%u(k)%inv = 0.0
147                      iv%tamdar(n)%v(k)%inv = 0.0
148                   endif
149                endif
150             endif
152          failed=.false.
153          if( iv%tamdar(n)%t(k)%qc >= obs_qc_pointer )  then
154          call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%t(k), max_error_t ,failed)
155          if( iv%info(tamdar)%proc_domain(k,n) ) then
156              num_qcstat_conv(1,tamdar,3,ipr) = num_qcstat_conv(1,tamdar,3,ipr) + 1
157          if(failed) then
158           num_qcstat_conv(2,tamdar,3,ipr) = num_qcstat_conv(2,tamdar,3,ipr) + 1
159            if ( write_rej_obs_conv ) then
160            write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
161            'tamdar',ob_vars(3),iv%info(tamdar)%lat(k,n),iv%info(tamdar)%lon(k,n),0.01*iv%tamdar(n)%p(k)
162            end if
163          end if
164          end if
165          end if
167          failed=.false.
168          if( iv%tamdar(n)%q(k)%qc >= obs_qc_pointer ) then
169           if( iv%tamdar(n)%t(k)%qc == fails_error_max ) then
170           failed=.true.
171           iv%tamdar(n)%q(k)%qc  = fails_error_max
172           iv%tamdar(n)%q(k)%inv = 0.0
173           else
174           call da_max_error_qc (it,iv%info(tamdar), n, iv%tamdar(n)%q(k), max_error_q ,failed)
175           endif
176          if( iv%info(tamdar)%proc_domain(k,n) ) then
177                     num_qcstat_conv(1,tamdar,4,ipr) = num_qcstat_conv(1,tamdar,4,ipr) + 1
178          if(failed) then
179          num_qcstat_conv(2,tamdar,4,ipr) = num_qcstat_conv(2,tamdar,4,ipr) + 1
180            if ( write_rej_obs_conv ) then
181            write(qcstat_conv_unit,'(2x,a10,2x,a4,3f12.2)')&
182            'tamdar',ob_vars(4),iv%info(tamdar)%lat(k,n),iv%info(tamdar)%lon(k,n),0.01*iv%tamdar(n)%p(k)
183            end if
184          end if
185          end if
186          end if
188       end do
189    end do
191    if (trace_use_dull) call da_trace_exit("da_check_max_iv_tamdar")
193 end subroutine da_check_max_iv_tamdar