da_transform_xtoy_metar_adj.inc

References to this file elsewhere.
1 subroutine da_transform_xtoy_metar_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          ! first guess state.
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 
16    real, allocatable :: model_u(:,:)
17    real, allocatable :: model_v(:,:)
18    real, allocatable :: model_t(:,:)
19    real, allocatable :: model_q(:,:)
20    real, allocatable :: model_psfc(:)
21 
22    if (trace_use_dull) call da_trace_entry("da_transform_xtoy_metar_adj")
23 
24    if (sfc_assi_options == sfc_assi_options_1) then
25       allocate (model_u(1,iv%info(metar)%n1:iv%info(metar)%n2))
26       allocate (model_v(1,iv%info(metar)%n1:iv%info(metar)%n2))
27       allocate (model_t(1,iv%info(metar)%n1:iv%info(metar)%n2))
28       allocate (model_q(1,iv%info(metar)%n1:iv%info(metar)%n2))
29       allocate (model_psfc(iv%info(metar)%n1:iv%info(metar)%n2))
30       ! [1.2] Interpolate horizontally:
31       do n=iv%info(metar)%n1, iv%info(metar)%n2
32          model_u(1,n)  = jo_grad_y%metar(n)%u
33          model_v(1,n)  = jo_grad_y%metar(n)%v
34          model_t(1,n)  = jo_grad_y%metar(n)%t
35          model_q(1,n)  = jo_grad_y%metar(n)%q
36          model_psfc(n) = jo_grad_y%metar(n)%p
37       end do
38       call da_interp_lin_3d_adj (jo_grad_x%u, iv%info(metar), model_u)
39       call da_interp_lin_3d_adj (jo_grad_x%v, iv%info(metar), model_v)
40       call da_interp_lin_3d_adj (jo_grad_x%t, iv%info(metar), model_t)
41       call da_interp_lin_3d_adj (jo_grad_x%q, iv%info(metar), model_q)
42       call da_interp_lin_2d_adj (jo_grad_x%psfc, iv%info(metar), 1, model_psfc)
43 
44       deallocate (model_u)
45       deallocate (model_v)
46       deallocate (model_t)
47       deallocate (model_q)
48       deallocate (model_psfc)
49    else if (sfc_assi_options == sfc_assi_options_2) then
50       call da_transform_xtopsfc_adj(grid,iv,metar,iv%metar(:), jo_grad_y%metar(:),jo_grad_x)
51    end if
52 
53    if (trace_use_dull) call da_trace_exit("da_transform_xtoy_metar_adj")
54 
55 end subroutine da_transform_xtoy_metar_adj
56 
57