wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_varbc / da_varbc_tl.inc
blobeef1a10ce217c2a8fb1cb77615445fc43b98c63f
1   subroutine da_varbc_tl (cv_size, cv, iv, y)
3    !---------------------------------------------------------------------------
4    !  PURPOSE: Apply bias correction to radiance innovations.
5    !
6    !  METHOD: delta_d    = y - delta_bc
7    !          y = H (delta_x)
8    !          delta_bc   = SUM( delta_beta_i Pred_i )
9    !                       i
10    !          delta_beta = VarBC parameters
11    !          Pred = Bias predictors
12    ! 
13    !  Called from da_transform_vtoy
14    !
15    !  HISTORY: 10/27/2007 - Creation                     Tom Auligne
16    !---------------------------------------------------------------------------
18    implicit none
20    integer, intent(in)           :: cv_size         ! Size of cv array.
21    real, intent(in)              :: cv(1:cv_size)   ! control variables.
22    type (iv_type), intent(in)    :: iv              ! O-B structure.
23    type (y_type), intent(inout)  :: y               ! y = h (xa)
25    integer              :: inst, i, k, num_rad, n, npred
26    real, allocatable    :: varbc_param_tl(:)
27    
28    if (trace_use) call da_trace_entry("da_varbc_tl")
30       do inst = 1, iv%num_inst                 ! loop for sensor
31          num_rad = iv%instid(inst)%num_rad
32          if (num_rad < 1) cycle
34          allocate(varbc_param_tl(iv%instid(inst)%varbc_info%npredmax))
35   
36          do k=1,iv%instid(inst)%nchan          ! loop for channel
37             npred = iv%instid(inst)%varbc(k)%npred
38             if (npred <= 0) cycle              ! VarBC channels only
39                      
40            !---------------------------------------------------------------
41            ! Change of variable (preconditioning) 
42            !---------------------------------------------------------------
43             varbc_param_tl(:) = 0.0
44             do i = 1, npred
45                varbc_param_tl(i) = SUM( cv(iv%instid(inst)%varbc(k)%index(1:npred)) * &
46                                         iv%instid(inst)%varbc(k)%vtox(i,1:npred) )
47             end do      
48                 
49            !---------------------------------------------------------------
50            ! TL of bias correction through linear regression
51            !---------------------------------------------------------------
52             do n = 1, num_rad                                         ! loop for pixel      
53                if (iv%instid(inst)%tb_qc(k,n) <= qc_varbc_bad) cycle  ! good obs only
55                 y%instid(inst)%tb(k,n) = y%instid(inst)%tb(k,n) +   &
56                    SUM( varbc_param_tl(1:npred) *                   &
57                         iv%instid(inst)%varbc_info%pred(            &
58                         iv%instid(inst)%varbc(k)%ipred(1:npred),n) )
59             end do
60          end do
62          deallocate(varbc_param_tl)
63       end do                                   ! end loop for sensor
65    if (trace_use) call da_trace_exit("da_varbc_tl")
67  end subroutine da_varbc_tl