updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_varbc_tamdar / da_varbc_tamdar_adj.inc
blob0ebf7d5353b15829241fba1587315be697915220
1    subroutine da_varbc_tamdar_adj (cv_size, cv, iv, y)
3    implicit none
5    integer, intent(in)           :: cv_size    
6    real, intent(inout)           :: cv(1:cv_size)
7    type (iv_type), intent(inout) :: iv        
8    type (y_type), intent(in)     :: y             ! y = h (xa)
10    integer                       :: isn,iflt,iobs,ipred,iphase,npred
11    real, allocatable             :: varbc_param_adj(:)
12    real                          :: yh, cv_local
15    if (trace_use) call da_trace_entry("da_varbc_tamdar_adj")
17    npred = iv%varbc_tamdar%npred
19    allocate( varbc_param_adj(npred) )
21    do iflt = 1, iv%varbc_tamdar%nair
22       do iphase = 1, iv%varbc_tamdar%nphase
23          if (iv%varbc_tamdar%nobs_sum(iphase,iflt) >= varbc_tamdar_nobsmin) then
25              varbc_param_adj(:) = 0.0
27              if (iv%varbc_tamdar%nobs(iphase,iflt) > 0 .and. iv%varbc_tamdar%ifuse(iphase,iflt) > 0) then
28                  do iobs = 1, iv%varbc_tamdar%nobs(iphase,iflt)
29                     isn = iv%varbc_tamdar%obs_sn(iobs,iphase,iflt)
30                     yh = 0
31                     if( iv%tamdar(isn)%t(1)%qc >= 0 ) yh = y%tamdar(isn)%t(1)
32                     varbc_param_adj(1:npred) = varbc_param_adj(1:npred) + &
33                                                yh * iv%varbc_tamdar%pred(1:npred,iphase,iflt)
34                  end do
35              end if
37              do ipred = 1, npred
38                 cv_local = SUM( varbc_param_adj(1:npred) * iv%varbc_tamdar%vtox(ipred,1:npred,iphase,iflt) )
39                 cv_local = wrf_dm_sum_real(cv_local)
41                 if (iv%varbc_tamdar%nobs(iphase,iflt) > 0 .and. iv%varbc_tamdar%ifuse(iphase,iflt) > 0) &
42                     cv(iv%varbc_tamdar%index(ipred,iphase,iflt)) = cv( iv%varbc_tamdar%index(ipred,iphase,iflt) ) + cv_local
43              end do
45          end if
46       end do
47    end do
49    deallocate(varbc_param_adj)       
51    if (trace_use) call da_trace_exit("da_varbc_tamdar_adj")
53    end subroutine da_varbc_tamdar_adj