da_transfer_wrftltoxa.inc

References to this file elsewhere.
1 subroutine da_transfer_wrftltoxa(grid, config_flags, filnam )
2 
3    !---------------------------------------------------------------------------
4    ! Purpose: Convert WRFTL variables to analysis increments
5    !           (inverse of the incremental part of xatowrf)
6    !---------------------------------------------------------------------------
7 
8    implicit none
9    
10    type(domain),               intent(inout) :: grid
11    type(grid_config_rec_type), intent(inout) :: config_flags
12 
13    character*4,                intent(in)    :: filnam
14 
15    integer :: i, j, k
16    integer :: is, ie, js, je, ks, ke
17 
18    real :: sdmd, s1md
19    real :: g_press(grid%xp%ims:grid%xp%ime,grid%xp%jms:grid%xp%jme, &
20       grid%xp%kms:grid%xp%kme)
21 
22    integer ndynopt
23 
24    if (trace_use) call da_trace_entry("da_transfer_wrftltoxa")
25 
26    !---------------------------------------------------------------------------
27    ! [0.0] input
28    !---------------------------------------------------------------------------
29 
30    is=grid%xp%its
31    ie=grid%xp%ite
32    js=grid%xp%jts
33    je=grid%xp%jte
34    ks=grid%xp%kts
35    ke=grid%xp%kte
36 
37    ndynopt      = grid%dyn_opt
38    grid%dyn_opt = DYN_EM_TL
39    call nl_set_dyn_opt (1 , DYN_EM_TL)
40 
41    call da_med_initialdata_input(grid , config_flags, filnam)
42 
43    grid%dyn_opt = ndynopt
44    call nl_set_dyn_opt (1 , DYN_EM)
45 
46    !---------------------------------------------------------------------------
47    ! [1.0] Get the specific humidity increments from mixing ratio increments
48    !---------------------------------------------------------------------------
49 
50    do k=ks,ke
51       do j=js,je
52          do i=is,ie
53             grid%xa%q(i,j,k)=grid%g_moist(i,j,k,P_G_QV)*(1.0-grid%xb%q(i,j,k))**2
54          end do
55       end do
56    end do
57 
58    !---------------------------------------------------------------------------
59    ! [2.0] COMPUTE psfc increments from mu-increments
60    !---------------------------------------------------------------------------
61 
62    do j=js,je
63       do i=is,ie
64          sdmd=0.0
65          s1md=0.0
66          do k=ks,ke
67             sdmd=sdmd+grid%g_moist(i,j,k,P_G_QV)*grid%em_dnw(k)
68             s1md=s1md+(1.0+grid%moist(i,j,k,P_QV))*grid%em_dnw(k)
69          end do
70          grid%xa%psfc(i,j)=-grid%xb%psac(i,j)*sdmd-grid%em_g_mu_2(i,j)*s1md
71       end do
72    end do
73 
74    !---------------------------------------------------------------------------
75    ! [3.0] COMPUTE pressure increments 
76    !---------------------------------------------------------------------------
77 
78    do j=js,je
79       do i=is,ie
80          g_press(i,j,ke+1)=0.0
81          do k=ke,ks,-1
82             g_press(i,j,k)=g_press(i,j,k+1) &
83                -(grid%em_g_mu_2(i,j)*(1.0+grid%moist(i,j,k,P_QV)) &
84                +(grid%em_mu_2(i,j)+grid%em_mub(i,j))* &
85                grid%g_moist(i,j,k,P_G_QV))*grid%em_dn(k)
86             grid%xa%p(i,j,k)=0.5*(g_press(i,j,k)+g_press(i,j,k+1))
87          end do
88       end do
89    end do
90 
91    !---------------------------------------------------------------------------
92    ! [4.0] convert theta increments to t increments
93    !---------------------------------------------------------------------------
94 
95    do k=ks,ke
96       do j=js,je
97          do i=is,ie
98             grid%xa%t(i,j,k)=grid%xb%t(i,j,k)*(grid%em_g_t_2(i,j,k)/ &
99                (300.0+grid%em_t_2(i,j,k)) &
100                +kappa*grid%xa%p(i,j,k)/grid%xb%p(i,j,k))
101          end do
102       end do
103    end do
104 
105    ! FIX. In the inverse, g_ph information is lost. This should be investigated 
106    ! later.
107 
108    !-------------------------------------------------------------------------
109    ! [5.0] convert from c-grid to a-grid
110    !-------------------------------------------------------------------------
111 
112    ! fill the halo region for u and v.
113    call wrf_dm_halo(grid%xp%domdesc,grid%xp%comms,grid%xp%halo_id3)
114 
115    do k=ks,ke
116       do j=js,je
117          do i=is,ie
118             grid%xa%u(i,j,k)=0.5*(grid%em_g_u_2(i+1,j,  k)+grid%em_g_u_2(i,j,k))
119             grid%xa%v(i,j,k)=0.5*(grid%em_g_v_2(i  ,j+1,k)+grid%em_g_v_2(i,j,k))
120          end do
121       end do
122    end do
123 
124    !---------------------------------------------------------------------------
125    ! [6.0] all the simple ones
126    !---------------------------------------------------------------------------
127 
128    do j=js,je
129       do k=ks,ke+1
130          do i=is,ie
131             grid%xa%w(i,j,k)=grid%em_g_w_2(i,j,k)
132          end do
133       end do
134    end do
135 
136 #ifdef VAR4D_MICROPHYSICS
137    ! New code needed once we introduce the microphysics to 4dvar in 2008
138    if (size(grid%moist,dim=4) >= 4) then
139       do k=ks,ke
140          do j=js,je
141             do i=is,ie
142                grid%xa%qcw(i,j,k)=grid%g_moist(i,j,k,p_qcw)
143                grid%xa%qrn(i,j,k)=grid%g_moist(i,j,k,p_qrn)
144             end do
145          end do
146       end do
147    end if
148 
149    if (size(grid%moist,dim=4) >= 6) then
150       do k=ks,ke
151          do j=js,je
152             do i=is,ie
153                grid%xa%qci(i,j,k)=grid%g_moist(i,j,k,p_qci)
154                grid%xa%qsn(i,j,k)=grid%g_moist(i,j,k,p_qsn)
155             end do
156          end do
157       end do
158    end if
159 
160    if (size(grid%moist_2,dim=4) >= 7) then
161       do k=ks,ke
162          do j=js,je
163             do i=is,ie
164                grid%xa%qgr(i,j,k)=grid%g_moist(i,j,k,p_qgr) 
165             end do
166          end do
167       end do
168    end if
169 #endif
170 
171    if (trace_use) call da_trace_exit("da_transfer_wrftltoxa")
172 
173 end subroutine da_transfer_wrftltoxa
174 
175