Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_varbc / da_varbc_direct.inc
blob8c4ca76c6790f12d6501c091b443adce0c4a110d
1   subroutine da_varbc_direct (iv)
3    !---------------------------------------------------------------------------
4    !  PURPOSE: Apply bias correction to radiance innovations.
5    !
6    !  METHOD:  d   = y - H(x) - bc
7    !          bc   = SUM( beta_i Pred_i )
8    !                  i
9    !          beta = VarBC parameters
10    !          Pred = Bias predictors
11    !
12    !  Called from da_get_innov_vector_radiance
13    !
14    !  HISTORY: 10/26/2007 - Creation                     Tom Auligne
15    !---------------------------------------------------------------------------
17    implicit none
19    type (iv_type), intent(inout)  :: iv        ! O-B structure.
21    integer              :: inst, k, n, npred, npredmax
22     
23    if (trace_use) call da_trace_entry("da_varbc_direct")
25       do inst = 1, iv%num_inst                 ! loop for sensor
26          npredmax = iv%instid(inst)%varbc_info%npredmax
27          write(unit=stdout,fmt='(A,A)') 'VARBC: Applying bias correction for ', &
28             trim(iv%instid(inst)%rttovid_string)
29                   
30          do k=1,iv%instid(inst)%nchan          ! loop for channel
31             npred          = iv%instid(inst)%varbc(k)%npred
32             if (npred <= 0) cycle                      ! VarBC channel only
34             do n= iv%instid(inst)%info%n1,iv%instid(inst)%info%n2  ! loop for pixel                  
35               ! apply bias correction through linear regression
36               !-------------------------------------------------
37               if ( iv%instid(inst)%tb_inv(k,n) > missing_r ) then
38                  iv%instid(inst)%tb_inv(k,n) = iv%instid(inst)%tb_inv(k,n) - &
39                        SUM( iv%instid(inst)%varbc(k)%param(1:npred) * &
40                             iv%instid(inst)%varbc_info%pred( &
41                             iv%instid(inst)%varbc(k)%ipred(1:npred),n) )
42               end if
43             end do
44          end do
45       end do                                   ! end loop for sensor
47    if (trace_use) call da_trace_exit("da_varbc_direct")
49  end subroutine da_varbc_direct