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