da_transform_xtoy_satem_adj.inc

References to this file elsewhere.
1 subroutine da_transform_xtoy_satem_adj(iv, xp, jo_grad_y, jo_grad_x, xb)
2 
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
6 
7    implicit none
8 
9    type (ob_type), intent(in)    :: iv          ! obs. inc vector (o-b).
10    type (xpose_type), intent(in) :: xp          ! Domain decomposition vars.
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    type (xb_type),  intent(in)   :: xb          ! first guess state.
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)       :: ADJ_pre_ma,ADJ_tv_ma
27    real                           :: ADJ_tv1,ADJ_tv2
28 
29    ADJ_pre_ma(:) = 0.
30    ADJ_tv_ma(:)  = 0.
31    ADJ_tv1 = 0.
32    ADJ_tv2 = 0.
33 
34    if (iv%num_satem > 0) then
35       do n=iv%ob_numb(iv%current_ob_time-1)%satem + 1, &
36          iv%ob_numb(iv%current_ob_time)%satem
37          num_levs = iv % satem(n) % info % levels
38 
39          ! [1.0] Get horizontal interpolation weights:
40 
41          i = iv%satem(n)%loc%i
42          dy = iv%satem(n)%loc%dy
43          dym = iv%satem(n)%loc%dym
44          j = iv%satem(n)%loc%j
45          dx = iv%satem(n)%loc%dx
46          dxm = iv%satem(n)%loc%dxm
47          ks = xp%kts; ke = xp%kte
48 
49          ! [2.0] Virtual temperature at obs pt.
50 
51          call da_tv_profile(xp,xb,i,j,dx,dxm,dy,dym,pre_ma,tv_ma)
52 
53          ! [3.0] Find model vertical position of pressure and do interp.
54 
55          call da_find_layer(layer2,tv2,iv%satem(n)%ref_p,pre_ma,tv_ma,ks,ke)
56          pres2 = iv%satem(n)%ref_p
57 
58          ! [4.0] Adjoint calculation of Satem thickness
59 
60          do k=1, num_levs
61             if (ABS(iv % satem(n) %p (k) - missing_r) > 1.) then
62                call da_find_layer(layer1,tv1,iv%satem(n)%p(k),pre_ma,tv_ma, &
63                  ks,ke)
64 
65                call da_thickness_adj(pre_ma,tv_ma,ks,ke,tv1,tv2,layer1, &
66                    layer2, iv%satem(n)%p(k),pres2,ADJ_pre_ma,ADJ_tv_ma, &
67                    ADJ_tv1,ADJ_tv2,jo_grad_y%satem(n)%thickness(k))
68 
69                call da_find_layer_adj(layer1,tv1,iv%satem(n)%p(k),         &
70                    pre_ma,tv_ma,ks,ke,ADJ_tv1,ADJ_pre_ma,ADJ_tv_ma)
71 
72                pres2 = iv%satem(n)%p(k)
73                layer2 = layer1
74                tv2 = tv1
75             end if
76          end do
77 
78          ! [5.0] Adjoint of layer-finding and vertical interpolation
79 
80          call da_find_layer_adj(layer2,tv2,iv%satem(n)%ref_p,              &
81               pre_ma,tv_ma,ks,ke,ADJ_tv2,ADJ_pre_ma,ADJ_tv_ma)
82 
83          ! [6.0] Adjoint of horizontal interpolation
84 
85          call da_tv_profile_adj(xp,jo_grad_x,xb,i,j,dx,dxm,dy,dym,         &
86               pre_ma,tv_ma,ADJ_pre_ma,ADJ_tv_ma)
87       end do
88    end if
89 
90 end subroutine da_transform_xtoy_satem_adj
91 
92