Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_physics / da_transform_xtopsfc.inc
blob7b938ac26e8fff0c15f9113d19871a79778200ba
1 subroutine da_transform_xtopsfc(grid, iv, obs_index, synop, y_synop)
3    !---------------------------------------------------------------------
4    ! Purpose: TBD
5    !---------------------------------------------------------------------
7    implicit none
9    type (domain),              intent(in)  :: grid
10    type (iv_type),             intent(in)  :: iv
11    integer,                    intent(in)  :: obs_index
12    type (synop_type),          intent(in)  :: synop(:)
13    type (residual_synop_type), intent(out) :: y_synop(:)
15    integer :: n
16    real :: to, qo
17    real, allocatable :: hsm(:,:)
18    real, allocatable :: tsm(:,:)
19    real, allocatable :: qsm(:,:)
20    real, allocatable :: psm(:,:)
21    real, allocatable :: psm_prime(:,:)
22    real, allocatable :: u(:,:)
23    real, allocatable :: v(:,:)
24    real, allocatable :: t(:,:)
25    real, allocatable :: q(:,:)
27    if (trace_use) call da_trace_entry("da_transform_xtopsfc")
29    allocate (hsm(1,iv%info(obs_index)%n1:iv%info(obs_index)%n2))
30    allocate (tsm(1,iv%info(obs_index)%n1:iv%info(obs_index)%n2))
31    allocate (qsm(1,iv%info(obs_index)%n1:iv%info(obs_index)%n2))
32    allocate (psm(1,iv%info(obs_index)%n1:iv%info(obs_index)%n2))
33    allocate (psm_prime(1,iv%info(obs_index)%n1:iv%info(obs_index)%n2))
34    allocate (u(1,iv%info(obs_index)%n1:iv%info(obs_index)%n2))
35    allocate (v(1,iv%info(obs_index)%n1:iv%info(obs_index)%n2))
36    allocate (t(1,iv%info(obs_index)%n1:iv%info(obs_index)%n2))
37    allocate (q(1,iv%info(obs_index)%n1:iv%info(obs_index)%n2))
39    ! 2.0 Surface assimilation approach 2 (2-m t and q, 10-m u and v)
41    call da_interp_lin_2d(grid%xa%u10,  iv%info(obs_index), 1, u)
42    call da_interp_lin_2d(grid%xa%v10,  iv%info(obs_index), 1, v)
43    call da_interp_lin_2d(grid%xa%psfc, iv%info(obs_index), 1, psm_prime)
44    call da_interp_lin_2d(grid%xa%t2,   iv%info(obs_index), 1, t)
45    call da_interp_lin_2d(grid%xa%q2,   iv%info(obs_index), 1, q)
46    do n=iv%info(obs_index)%n1,iv%info(obs_index)%n2
47       y_synop(n)%u=u(1,n)
48       y_synop(n)%v=v(1,n)
49       y_synop(n)%t=t(1,n)
50       y_synop(n)%q=q(1,n)
51    end do
53    ! 3.2 model background surface p, t, q, h at observed site:
55    call da_interp_lin_2d(grid%xb%terr, iv%info(obs_index), 1, hsm)
56    call da_interp_lin_2d(grid%xb%t2,   iv%info(obs_index), 1, tsm)
57    call da_interp_lin_2d(grid%xb%q2,   iv%info(obs_index), 1, qsm)
58    call da_interp_lin_2d(grid%xb%psfc, iv%info(obs_index), 1, psm)
60    do n=iv%info(obs_index)%n1,iv%info(obs_index)%n2
61       if (synop(n)%p%qc < 0) then
62          y_synop(n)%p = 0.0
63       else
65          ! 3.0 The pressure at the observed height: 
67          ! 3.1 Constants:
69           to = -888888.0
70           qo = -888888.0
71              
72          ! Terrain height at the observed site:
74          ! 3.3 perturbations t, q, p at the model surface:
76          ! 4.0 Compute the perturbation of the surface pressure perturbation 
77          !     at the observed height
79          if (synop(n)%t%qc >= 0 .and. synop(n)%q%qc >= 0) then
80             ! 4.1 Observed value = background + innovation: both t, q available:
81             !     ----------------------------------------
83             to = tsm(1,n) + synop(n)%t%inv
84             qo = qsm(1,n) + synop(n)%q%inv
85             call da_sfc_pre_lin(y_synop(n)%p, psm_prime(1,n), y_synop(n)%t, y_synop(n)%q, &
86                psm(1,n), tsm(1,n), qsm(1,n), hsm(1,n), synop(n)%h, to, qo)
87          else if (synop(n)%t%qc >= 0 .and. synop(n)%q%qc < 0) then
89             ! 4.2 Observed value = background + innovation: only t available
90             !     ----------------------------------------
92             to = tsm(1,n) + synop(n)%t%inv
93             call da_sfc_pre_lin(y_synop(n)%p, psm_prime(1,n), y_synop(n)%t, y_synop(n)%q, &
94                psm(1,n), tsm(1,n), qsm(1,n), hsm(1,n), synop(n)%h, to)
95          else
96             ! 4.3 No observed t and q available:
97             !     -----------------------------
99             call da_sfc_pre_lin(y_synop(n)%p, psm_prime(1,n), y_synop(n)%t, y_synop(n)%q, &
100                psm(1,n), tsm(1,n), qsm(1,n), hsm(1,n), synop(n)%h)
101          end if
102       end if
103    end do
105    deallocate (hsm)
106    deallocate (tsm)
107    deallocate (qsm)
108    deallocate (psm)
109    deallocate (psm_prime)
110    deallocate (u)
111    deallocate (v)
112    deallocate (t)
113    deallocate (q)
115    if (trace_use) call da_trace_exit("da_transform_xtopsfc")
117 end subroutine da_transform_xtopsfc