da_tv_profile_tl.inc
References to this file elsewhere.
1 subroutine da_tv_profile_tl(grid,i,j,dx,dxm,dy,dym,pre_ma,tv_ma,TGL_pre_ma,TGL_tv_ma)
2
3 !--------------------------------------------------------------------------
4 ! Purpose: tangent-linear routine for da_tv_profile
5 !--------------------------------------------------------------------------
6
7 implicit none
8
9 type (domain), intent(in) :: grid
10 integer, intent(in) :: i, j ! OBS location
11 real, intent(in) :: dx, dxm ! interpolation weights.
12 real, intent(in) :: dy, dym ! interpolation weights.
13 real, intent(out) :: TGL_pre_ma(kts-1:kte+1)
14 real, intent(out) :: TGL_tv_ma(kts-1:kte+1)
15 real, intent(out) :: pre_ma(kts-1:kte+1)
16 real, intent(out) :: tv_ma(kts-1:kte+1)
17
18 integer :: ii,jj
19 real :: tv_m(2,2,kts:kte)
20 real :: TGL_tv_m(2,2,kts:kte)
21
22 if (trace_use_dull) call da_trace_entry("da_tv_profile_tl")
23
24 do ii=i,i+1
25 do jj=j,j+1
26 TGL_tv_m(ii-i+1,jj-j+1,kts:kte) = grid%xa%t(ii,jj,kts:kte)*(1.0+0.61*grid%xb%q(ii,jj,kts:kte)) + &
27 0.61*grid%xb%t(ii,jj,kts:kte)*grid%xa%q(ii,jj,kts:kte)
28 tv_m(ii-i+1,jj-j+1,kts:kte) = grid%xb%t(ii,jj,kts:kte)*(1.0+0.61*grid%xb%q(ii,jj,kts:kte))
29 end do
30 end do
31
32 TGL_pre_ma(kts:kte) = dym* ( dxm * grid%xa%p(i,j,kts:kte) + dx * grid%xa%p(i+1,j,kts:kte) ) + &
33 dy * ( dxm * grid%xa%p(i,j+1,kts:kte) + dx * grid%xa%p(i+1,j+1,kts:kte) )
34 TGL_tv_ma (kts:kte) = dym* ( dxm * TGL_tv_m(1,1,kts:kte) + dx * TGL_tv_m(2,1,kts:kte) ) + &
35 dy * ( dxm * TGL_tv_m(1,2,kts:kte) + dx * TGL_tv_m(2,2,kts:kte) )
36 pre_ma(kts:kte) = dym* ( dxm * grid%xb%p(i,j,kts:kte) + dx * grid%xb%p(i+1,j,kts:kte) ) + &
37 dy * ( dxm * grid%xb%p(i,j+1,kts:kte) + dx * grid%xb%p(i+1,j+1,kts:kte) )
38 tv_ma (kts:kte) = dym* ( dxm * tv_m (1,1,kts:kte) + dx * tv_m (2,1,kts:kte) ) + &
39 dy * ( dxm * tv_m (1,2,kts:kte) + dx * tv_m (2,2,kts:kte) )
40
41 if (trace_use_dull) call da_trace_exit("da_tv_profile_tl")
42
43 end subroutine da_tv_profile_tl
44
45