da_transform_xtoy_satem_adj.inc
References to this file elsewhere.
1 subroutine da_transform_xtoy_satem_adj(grid, iv, jo_grad_y, jo_grad_x)
2
3 !-----------------------------------------------------------------------
4 ! Purpose: TBD
5 !-----------------------------------------------------------------------
6
7 implicit none
8
9 type (domain), intent(in) :: grid
10 type (iv_type), intent(in) :: iv ! obs. inc vector (o-b).
11 type (y_type) , intent(inout) :: jo_grad_y ! grad_y(jo)
12 type (x_type) , intent(inout) :: jo_grad_x ! grad_x(jo)
13
14 integer :: n ! Loop counter.
15 integer :: num_levs ! obs vertical levels
16
17 integer :: k
18 real :: pre_ma(kts-1:kte+1)
19 real :: tv_ma(kts-1:kte+1)
20 integer :: layer1,layer2
21 real :: tv1,tv2,pres2
22
23 real :: ADJ_pre_ma(kts-1:kte+1)
24 real :: ADJ_tv_ma(kts-1:kte+1)
25 real :: ADJ_tv1,ADJ_tv2
26
27 if (trace_use_dull) call da_trace_entry("da_transform_xtoy_satem_adj")
28
29 ADJ_pre_ma(:) = 0.0
30 ADJ_tv_ma(:) = 0.0
31 ADJ_tv1 = 0.0
32 ADJ_tv2 = 0.0
33
34 do n=iv%info(satem)%n1,iv%info(satem)%n2
35 num_levs = iv%info(satem)%levels(n)
36
37 ! [2.0] Virtual temperature at obs pt.
38
39 call da_tv_profile(grid,iv%info(satem),n,pre_ma,tv_ma)
40
41 ! [3.0] Find model vertical position of pressure and do interp.
42
43 call da_find_layer(layer2,tv2,iv%satem(n)%ref_p,pre_ma,tv_ma,kts,kte)
44 pres2 = iv%satem(n)%ref_p
45
46 ! [4.0] Adjoint calculation of satem thickness
47
48 do k=1, num_levs
49 if (ABS(iv % satem(n) %p (k) - missing_r) > 1.0) then
50 call da_find_layer(layer1,tv1,iv%satem(n)%p(k),pre_ma,tv_ma, &
51 kts,kte)
52
53 call da_thickness_adj(pre_ma,tv_ma,kts,kte,tv1,tv2,layer1, &
54 layer2, iv%satem(n)%p(k),pres2,ADJ_pre_ma,ADJ_tv_ma, &
55 ADJ_tv1,ADJ_tv2,jo_grad_y%satem(n)%thickness(k))
56
57 call da_find_layer_adj(layer1,tv1,iv%satem(n)%p(k), &
58 pre_ma,tv_ma,kts,kte,ADJ_tv1,ADJ_pre_ma,ADJ_tv_ma)
59
60 pres2 = iv%satem(n)%p(k)
61 layer2 = layer1
62 tv2 = tv1
63 end if
64 end do
65
66 ! [5.0] Adjoint of layer-finding and vertical interpolation
67
68 call da_find_layer_adj(layer2,tv2,iv%satem(n)%ref_p, &
69 pre_ma,tv_ma,kts,kte,ADJ_tv2,ADJ_pre_ma,ADJ_tv_ma)
70
71 ! [6.0] Adjoint of horizontal interpolation
72 call da_tv_profile_adj(grid,jo_grad_x,iv%info(satem),n, &
73 pre_ma,tv_ma,ADJ_pre_ma,ADJ_tv_ma)
74 end do
75
76 if (trace_use_dull) call da_trace_exit("da_transform_xtoy_satem_adj")
77
78 end subroutine da_transform_xtoy_satem_adj
79
80