da_transform_xtoy_satem.inc
References to this file elsewhere.
1 subroutine da_transform_xtoy_satem (grid, iv, y)
2
3 !-----------------------------------------------------------------------
4 ! Purpose: TBD
5 !-----------------------------------------------------------------------
6
7 implicit none
8
9 type (domain), intent(in) :: grid
10 type (iv_type), intent(in) :: iv ! Innovation vector (O-B).
11 type (y_type), intent(inout) :: y ! y = h (grid%xa)
12
13 integer :: n ! Loop counter.
14 integer :: i, j ! Index dimension.
15 real :: dx, dxm !
16 real :: dy, dym !
17 integer :: num_levs ! obs vertical levels
18
19 integer :: k
20 real :: pre_ma(kts-1:kte+1)
21 real :: tv_ma(kts-1:kte+1)
22 integer :: layer1,layer2
23 real :: tv1,tv2,pres2
24
25 real :: TGL_pre_ma(kts-1:kte+1)
26 real :: TGL_tv_ma(kts-1:kte+1)
27 real :: TGL_tv1,TGL_tv2
28
29 if (trace_use_dull) call da_trace_entry("da_transform_xtoy_satem")
30
31 do n=iv%info(satem)%n1,iv%info(satem)%n2
32
33 num_levs = iv%info(satem)%levels(n)
34
35 ! [1.0] Get horizontal interpolation weights:
36
37 i = iv%info(satem)%i(1,n)
38 dy = iv%info(satem)%dy(1,n)
39 dym = iv%info(satem)%dym(1,n)
40 j = iv%info(satem)%j(1,n)
41 dx = iv%info(satem)%dx(1,n)
42 dxm = iv%info(satem)%dxm(1,n)
43
44 ! [2.0] Virtual temperature at obs pt.
45
46 call da_tv_profile_tl(grid,i,j,dx,dxm,dy,dym, &
47 pre_ma,tv_ma,TGL_pre_ma,TGL_tv_ma)
48
49 ! [3.0] Find model vertical position of pressure and do interp.
50
51 call da_find_layer_tl(layer2,tv2,iv%satem(n)%ref_p, &
52 pre_ma,tv_ma,kts,kte,TGL_tv2,TGL_pre_ma,TGL_tv_ma)
53 pres2 = iv%satem(n)%ref_p
54
55 ! [4.0] Thickness calculation
56
57 do k=1, num_levs
58 if (ABS(iv % satem(n) %p (k) - missing_r) > 1.0) then
59 call da_find_layer_tl(layer1,tv1,iv%satem(n)%p(k), &
60 pre_ma,tv_ma,kts,kte,TGL_tv1,TGL_pre_ma,TGL_tv_ma)
61 call da_thickness_tl(pre_ma,tv_ma,kts,kte,tv1,tv2,layer1,layer2,&
62 iv%satem(n)%p(k),pres2,TGL_pre_ma,TGL_tv_ma, &
63 TGL_tv1,TGL_tv2,y%satem(n)%thickness(k))
64
65 pres2 = iv%satem(n)%p(k)
66 layer2 = layer1
67 tv2 = tv1
68 end if
69 end do
70 end do
71
72 if (trace_use_dull) call da_trace_exit("da_transform_xtoy_satem")
73
74 end subroutine da_transform_xtoy_satem
75
76