Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_varbc / da_varbc_adj.inc
blob6b758986d075051f27567ebc8f3a4692b9d7deb1
1   subroutine da_varbc_adj (cv_size, cv, iv, y)
3    !---------------------------------------------------------------------------
4    !  PURPOSE: Apply bias correction to radiance innovations.
5    !
6    !  Called from da_transform_vtoy_adj
7    !
8    !  HISTORY: 10/27/2007 - Creation                     Tom Auligne
9    !---------------------------------------------------------------------------
11    implicit none
13    integer, intent(in)           :: cv_size         ! Size of cv array.
14    real, intent(inout)           :: cv(1:cv_size)   ! control variables.
15    type (iv_type), intent(inout) :: iv              ! O-B structure.
16    type (y_type), intent(in)     :: y               ! y = h (xa)
18    integer              :: inst, i, k, num_rad, n, npred
19    real, allocatable    :: varbc_param_adj(:)
20    real                 :: cv_local
22    if (trace_use) call da_trace_entry("da_varbc_adj")
24       do inst = 1, iv%num_inst                 ! loop for sensor
25          num_rad = iv%instid(inst)%num_rad
26    
27          allocate(varbc_param_adj(iv%instid(inst)%varbc_info%npredmax))
29          do k=1,iv%instid(inst)%nchan          ! loop for channel
30             npred = iv%instid(inst)%varbc(k)%npred
31             if (npred <= 0) cycle              ! VarBC channels only
32            
33            !---------------------------------------------------------------
34            ! Adjoint to bias correction through linear regression
35            !---------------------------------------------------------------
36             varbc_param_adj(:) = 0.0
37             if (num_rad >0) then
38                do n = 1, num_rad                                     ! loop for pixel      
39                   if (iv%instid(inst)%tb_qc(k,n) <= qc_varbc_bad) cycle  ! good obs only
40                   if (.not. iv%instid(inst)%info%proc_domain(1,n))  cycle   ! do not sum up HALO data
41                
42                   varbc_param_adj(1:npred) = varbc_param_adj(1:npred)      + &
43                                              y%instid(inst)%tb(k,n)        * &
44                                             iv%instid(inst)%varbc_info%pred( &
45                                             iv%instid(inst)%varbc(k)%ipred(1:npred),n ) 
46                end do           
47             end if
48                        
49            !---------------------------------------------------------------
50            ! Change of variable (preconditioning) + sum across processors
51            !---------------------------------------------------------------
52             do i = 1, npred
53                cv_local = SUM(varbc_param_adj(1:npred) * iv%instid(inst)%varbc(k)%vtox(i,1:npred))
54                cv(iv%instid(inst)%varbc(k)%index(i)) = cv(iv%instid(inst)%varbc(k)%index(i)) + &
55                                                        wrf_dm_sum_real(cv_local)        
56             end do
57          end do
58          
59          deallocate(varbc_param_adj)         
60       end do                                   ! end loop for sensor
62    if (trace_use) call da_trace_exit("da_varbc_adj")
64  end subroutine da_varbc_adj