da_transform_xtoy_metar_adj.inc

References to this file elsewhere.
1 subroutine da_transform_xtoy_metar_adj(xb, iv, xp, jo_grad_y, jo_grad_x)
2 
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
6 
7    implicit none
8 
9    type (xb_type),    intent(in) :: xb          ! first guess state.
10    type (ob_type),    intent(in) :: iv          ! obs. inc vector (o-b).
11    type (xpose_type), intent(in) :: xp          ! Domain decomposition vars.
12    type (y_type) , intent(inout) :: jo_grad_y   ! grad_y(jo)
13    type (x_type) , intent(inout) :: jo_grad_x   ! grad_x(jo)
14 
15    integer                       :: n        ! Loop counter.
16    integer                       :: i, j     ! Index dimension.
17    real                          :: dx, dxm  ! 
18    real                          :: dy, dym  !
19 
20    if (trace_use) call da_trace_entry("da_transform_xtoy_metar_adj")
21 
22    if (iv%num_metar > 0) then
23       if (sfc_assi_options == 1) then
24          do n=iv%ob_numb(iv%current_ob_time-1)%metar + 1, &
25             iv%ob_numb(iv%current_ob_time)%metar
26 
27             ! [1.1] Get horizontal interpolation weights:
28 
29             i = iv%metar(n)%loc%i
30             dy = iv%metar(n)%loc%dy
31             dym = iv%metar(n)%loc%dym
32             j = iv%metar(n)%loc%j
33             dx = iv%metar(n)%loc%dx
34             dxm = iv%metar(n)%loc%dxm
35 
36             ! [1.2] Interpolate horizontally:
37             call da_interp_obs_lin_2d_adj(jo_grad_x % u, xp, i, j, dx, dy, &
38                dxm, dym, jo_grad_y%metar(n)%u, iv%metar(n)%zk)
39             call da_interp_obs_lin_2d_adj(jo_grad_x % v, xp, i, j, dx, dy, &
40                dxm, dym, jo_grad_y%metar(n)%v, iv%metar(n)%zk)
41             call da_interp_obs_lin_2d_adj(jo_grad_x % t, xp, i, j, dx, dy, &
42                dxm, dym, jo_grad_y%metar(n)%t, iv%metar(n)%zk)
43             call da_interp_obs_lin_2d_adj(jo_grad_x % q, xp, i, j, dx, dy, &
44                dxm, dym, jo_grad_y%metar(n)%q, iv%metar(n)%zk)
45 
46             call da_interp_lin_2d_adj(jo_grad_x % psfc, xp%ims, xp%ime, &
47                xp%jms, xp%jme, i, j, dx, dy, dxm, dym, jo_grad_y%metar(n)%p)
48          end do
49       else if (sfc_assi_options == 2) then
50          do n=iv%ob_numb(iv%current_ob_time-1)%metar + 1, &
51             iv%ob_numb(iv%current_ob_time)%metar
52             call da_transform_xtopsfc_adj(xb,xp,iv%metar(n), &
53                                             jo_grad_y%metar(n),jo_grad_x)
54          end do
55       end if
56    end if
57 
58    if (trace_use) call da_trace_exit("da_transform_xtoy_metar_adj")
59 
60 end subroutine da_transform_xtoy_metar_adj
61 
62