updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_varbc_tamdar / da_varbc_tamdar_tl.inc
blob9e843e860cdd6df9bed82e93810cb9e2659bd536
1    subroutine da_varbc_tamdar_tl (cv_size, cv, iv, y)
3    !---------------------------------------------------!
4    !  METHOD: delta_d    = y - delta_bc                !
5    !          y          = H (delta_x)                 !
6    !          delta_bc   = SUM( delta_beta_i Pred_i )  !
7    !---------------------------------------------------!
9    implicit none
11    integer, intent(in)           :: cv_size
12    real, intent(in)              :: cv(1:cv_size)
13    type (iv_type), intent(in)    :: iv   
14    type (y_type), intent(inout)  :: y               ! y = h (xa)
16    integer                       :: isn,ivar,iflt,ipred,iobs,iphase
17    integer                       :: npred
18    real                          :: delta_bc
20    real, allocatable             :: varbc_param_tl(:)
21    
23    if (trace_use) call da_trace_entry("da_varbc_tamdar_tl")
25    npred = iv%varbc_tamdar%npred
26    allocate( varbc_param_tl(npred) )
28    do iflt = 1, iv%varbc_tamdar%nair
29       do iphase = 1, iv%varbc_tamdar%nphase
31          if (iv%varbc_tamdar%nobs(iphase,iflt) > 0 .and. iv%varbc_tamdar%ifuse(iphase,iflt) > 0) then
33              varbc_param_tl(:) = 0.0
35              do ipred = 1, npred
36                 varbc_param_tl(ipred) = &
37                            SUM( cv(iv%varbc_tamdar%index(1:npred,iphase,iflt)) * &
38                                 iv%varbc_tamdar%vtox(ipred,1:npred,iphase,iflt) )
39              end do     
41              delta_bc = SUM( varbc_param_tl(1:npred)* &
42                              iv%varbc_tamdar%pred(1:npred,iphase,iflt) )
44              do iobs = 1, iv%varbc_tamdar%nobs(iphase,iflt)
45                 isn = iv%varbc_tamdar%obs_sn(iobs,iphase,iflt)
46                 if (iv%tamdar(isn)%t(1)%qc >= 0) &
47                     y%tamdar(isn)%t(1) = y%tamdar(isn)%t(1) + delta_bc
48              end do
50          end if
51       end do
52    end do
54    deallocate(varbc_param_tl)
56    if (trace_use) call da_trace_exit("da_varbc_tamdar_tl")
58    end subroutine da_varbc_tamdar_tl