da_transform_xtopsfc_adj.inc

References to this file elsewhere.
1 subroutine da_transform_xtopsfc_adj(xb, xp, synop, j_synop_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 (xpose_type), intent(in) :: xp       ! domain decomposition vars.
11 
12    type (synop_type),          intent(in)    :: synop
13    type (residual_synop_type), intent(inout) :: j_synop_y   ! grad_y(jo)
14    type (x_type),              intent(inout) :: jo_grad_x   ! grad_x(jo)
15    
16    integer                :: i, j     ! index dimension.
17    real                   :: dx, dxm  ! 
18    real                   :: dy, dym  !
19 
20    real    :: hsm, to, qo
21    real    :: tsm, qsm, psm
22    real    :: psm_prime
23 
24    !--------------------------------------------------------------------
25    ! 1.0 Get cross pt. horizontal interpolation weights:
26    !--------------------------------------------------------------------
27 
28    i   = synop%loc%i
29    dy  = synop%loc%dy
30    dym = synop%loc%dym
31    j   = synop%loc%j
32    dx  = synop%loc%dx
33    dxm = synop%loc%dxm
34 
35    psm_prime = 0.0
36    if (synop%p%qc >= 0) then
37 
38       !------------------------------------------------------------------------
39       ! 2.0 Calculate gradient with respect the pressure at the observed height: 
40       !------------------------------------------------------------------------
41 
42       ! 2.1 Terrain height at the observed site (xj, yi):
43 
44       call da_interp_lin_2d(xb%terr, xp%ims, xp%ime, xp%jms, xp%jme, &
45                          i, j, dx, dy, dxm, dym, hsm)
46       call da_interp_lin_2d(xb%t2, xp%ims, xp%ime, xp%jms, xp%jme, &
47                          i, j, dx, dy, dxm, dym, tsm)
48       call da_interp_lin_2d(xb%q2, xp%ims, xp%ime, xp%jms, xp%jme, &
49                          i, j, dx, dy, dxm, dym, qsm)
50       call da_interp_lin_2d(xb%psfc, xp%ims, xp%ime, xp%jms, xp%jme, &
51                          i, j, dx, dy, dxm, dym, psm)
52 
53       to = -888888.0
54       qo = -888888.0
55 
56       ! 2.3 Zero out the adjoint variables:
57 
58       !----------------------------------------------------------------
59       ! 3.0 Surface pressure gradient at the observed height
60       !----------------------------------------------------------------
61 
62       ! 4.0 Get the surface pressure gradient at the model surface height (hsm)
63       ! 4.1 Both observed to and qo available:
64       if (synop%t%qc >= 0 .and. synop%q%qc >= 0) then
65          to = tsm + synop%t%inv
66          qo = qsm + synop%q%inv
67          call da_sfc_pre_adj (j_synop_y%p, psm_prime, j_synop_y%t, j_synop_y%q, &
68                              psm, tsm, qsm, hsm, synop%h, to, qo)
69 
70          ! 4.2 only observed to available:
71       else if (synop%t%qc >= 0 .and. synop%q%qc < 0) then
72          to = tsm + synop%t%inv
73          call da_sfc_pre_adj (j_synop_y%p, psm_prime, j_synop_y%t, j_synop_y%q, &
74                              psm, tsm, qsm, hsm, synop%h, to)
75 
76          ! 4.3 Both observed to and qo missing:
77       else
78          call da_sfc_pre_adj (j_synop_y%p, psm_prime, j_synop_y%t, j_synop_y%q, &
79                              psm, tsm, qsm, hsm, synop%h)
80       end if
81 
82    end if
83 
84    ! 2.2 convert the jo_grad_y to the model grids (t2, q2, u10, v10, and psfc)
85    call da_interp_lin_2d_adj(jo_grad_x%t2, xp%ims, xp%ime, xp%jms, xp%jme, &
86                           i, j, dx, dy, dxm, dym, j_synop_y%t)
87    call da_interp_lin_2d_adj(jo_grad_x%q2, xp%ims, xp%ime, xp%jms, xp%jme, &
88                           i, j, dx, dy, dxm, dym, j_synop_y%q)
89    call da_interp_lin_2d_adj(jo_grad_x%u10, xp%ims, xp%ime, xp%jms, xp%jme, &
90                           i, j, dx, dy, dxm, dym, j_synop_y%u)
91    call da_interp_lin_2d_adj(jo_grad_x%v10, xp%ims, xp%ime, xp%jms, xp%jme, &
92                           i, j, dx, dy, dxm, dym, j_synop_y%v)
93    call da_interp_lin_2d_adj(jo_grad_x%psfc, xp%ims, xp%ime, xp%jms, xp%jme, &
94                           i, j, dx, dy, dxm, dym, psm_prime)
95 
96 end subroutine da_transform_xtopsfc_adj
97 
98