da_transfer_wrftltoxa_adj.inc

References to this file elsewhere.
1 subroutine da_transfer_wrftltoxa_adj(grid, config_flags, filnam)
2 
3    !---------------------------------------------------------------------------
4    ! Purpose: Convert analysis increments into WRFAD increments 
5    !          (following xatowrf, but only keep the increments)
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, ndynopt
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    real :: utmp(grid%xp%ims:grid%xp%ime,grid%xp%jms:grid%xp%jme, &
22       grid%xp%kms:grid%xp%kme)
23    real :: vtmp(grid%xp%ims:grid%xp%ime,grid%xp%jms:grid%xp%jme, &
24       grid%xp%kms:grid%xp%kme)
25 
26    if (trace_use) call da_trace_entry("da_transfer_wrftltoxa_adj")
27 
28    is=grid%xp%its
29    ie=grid%xp%ite
30    js=grid%xp%jts
31    je=grid%xp%jte
32    ks=grid%xp%kts
33    ke=grid%xp%kte
34 
35    grid%g_u_2  = 0.0
36    grid%g_v_2  = 0.0
37    grid%g_w_2  = 0.0
38    grid%g_t_2  = 0.0
39    grid%g_moist   = 0.0
40    grid%g_mu_2 = 0.0
41    grid%g_ph_2 = 0.0
42    !---------------------------------------------------------------------------
43    ! [6.0] ALL THE SIMPLE ONES
44    !---------------------------------------------------------------------------
45 
46 #ifdef DM_PARALLEL
47 #include "HALO_XA.inc"
48 #endif
49 
50    do k=ks,ke+1
51       do j=js,je
52          do i=is,ie
53             grid%g_w_2(i,j,k)=grid%g_w_2(i,j,k)+grid%xa%w(i,j,k)
54          end do
55       end do
56    end do
57    grid%xa%w(is:ie,js:je,ks:ke+1) = 0.0
58 
59 #ifdef VAR4D_MICROPHYSICS
60    ! New code needed once we introduce the microphysics to 4dvar in 2008
61    if (size(grid%moist,dim=4) >= 4) then
62       do k=ks,ke
63          do j=js,je
64             do i=is,ie
65                grid%g_moist(i,j,k,p_qcw) =  grid%xa%qcw(i,j,k)
66                grid%g_moist(i,j,k,p_qrn) =  grid%xa%qrn(i,j,k)
67             end do
68          end do
69       end do
70    end if
71 
72    if (size(grid%moist,dim=4) >= 6) then
73       do k=ks,ke
74          do j=js,je
75             do i=is,ie
76                grid%g_moist(i,j,k,p_qci) =  grid%xa%qci(i,j,k)
77                grid%g_moist(i,j,k,p_qsn) =  grid%xa%qsn(i,j,k)
78             end do
79          end do
80       end do
81    end if
82 
83    if (size(grid%moist,dim=4) >= 7) then
84       do k=ks,ke
85          do j=js,je
86             do i=is,ie
87                grid%g_moist(i,j,k,p_qgr) =  grid%xa%qgr(i,j,k)
88             end do
89          end do
90       end do
91    end if
92 
93 #endif
94 
95    !----------------------------------------------------------------------------
96    ! [5.0] convert from c-grid to a-grid
97    ! ----------------------------------------------------------------------------
98 
99 #ifdef DM_PARALLEL
100    utmp=grid%xa%u
101    vtmp=grid%xa%v
102 
103    ! The western boundary
104    if (is == grid%xp%ids) utmp(is-1,js:je,ks:ke)=0.0
105 
106    ! The southern boundary
107    if (js == grid%xp%jds) vtmp(is:ie,js-1,ks:ke)=0.0
108 
109    do k=ks,ke
110       do j=js,je
111          do i=is,ie
112             grid%g_u_2(i,j,k)=grid%g_u_2(i,j,k)+0.5*(utmp(i-1,j  ,k)+utmp(i,j,k))
113             grid%g_v_2(i,j,k)=grid%g_v_2(i,j,k)+0.5*(vtmp(i  ,j-1,k)+vtmp(i,j,k))
114          end do
115       end do
116    end do
117 
118    ! The eastern boundary
119    if (ie == grid%xp%ide)  &
120        grid%g_u_2(ie+1,js:je,ks:ke)=grid%g_u_2(ie+1,js:je,ks:ke)+grid%xa%u(ie,js:je,ks:ke)/2.0
121 
122    ! The northern boundary
123    if (je == grid%xp%jde)  &
124        grid%g_v_2(is:ie,je+1,ks:ke)=grid%g_v_2(is:ie,je+1,ks:ke)+grid%xa%v(is:ie,je,ks:ke)/2.0
125 
126 #else
127 
128    do k=ks,ke
129       do j=js,je
130          do i=is+1,ie
131             grid%g_u_2(i,j,k)=grid%g_u_2(i,j,k)+0.5*(grid%xa%u(i-1,j,k)+grid%xa%u(i,j,k))
132          end do
133       end do
134       do j=js+1,je
135          do i=is,ie
136             grid%g_v_2(i,j,k)=grid%g_v_2(i,j,k)+0.5*(grid%xa%v(i,j-1,k)+grid%xa%v(i,j,k))
137          end do
138       end do
139    end do
140 
141    ! The western boundary
142    grid%g_u_2(is,js:je,ks:ke)=grid%g_u_2(is,js:je,ks:ke)+grid%xa%u(is,js:je,ks:ke)/2.0
143 
144    ! The eastern boundary
145    grid%g_u_2(ie+1,js:je,ks:ke)=grid%g_u_2(ie+1,js:je,ks:ke)+grid%xa%u(ie,js:je,ks:ke)/2.0
146 
147    ! The southern boundary
148    grid%g_v_2(is:ie,js,ks:ke)=grid%g_v_2(is:ie,js,ks:ke)+grid%xa%v(is:ie,js,ks:ke)/2.0
149 
150    ! The northern boundary
151    grid%g_v_2(is:ie,je+1,ks:ke)=grid%g_v_2(is:ie,je+1,ks:ke)+grid%xa%v(is:ie,je,ks:ke)/2.0
152 
153 
154 #endif
155    grid%xa%u(is:ie, js:je, ks:ke) = 0.0
156    grid%xa%v(is:ie, js:je, ks:ke) = 0.0
157 
158    !---------------------------------------------------------------------------
159    ! [4.0] CONVERT THETA inCREMENTS TO T inCREMENTS
160    !---------------------------------------------------------------------------
161 
162    ! In the inverse, g_ph information is lost. This should be investigated later!!!
163    ! However, in this adjoint problem, a_ph should be set to 0.0 Otherwise, a_ph 
164    ! will be initialized randomly!
165 
166    grid%g_ph_2=0.0
167 
168    do k=ks,ke
169       do j=js,je
170          do i=is,ie
171             grid%xa%p(i,j,k)=grid%xa%p(i,j,k)+grid%xb%t(i,j,k)*kappa*grid%xa%t(i,j,k)/grid%xb%p(i,j,k)
172             grid%g_t_2(i,j,k)=grid%g_t_2(i,j,k)+grid%xb%t(i,j,k)*grid%xa%t(i,j,k)/(t0+grid%em_t_2(i,j,k))
173          end do
174       end do
175    end do
176 
177    grid%xa%t(is:ie, js:je, ks:ke) = 0.0
178    !---------------------------------------------------------------------------
179    ! [3.0] COMPUTE pressure increments 
180    !---------------------------------------------------------------------------
181 
182    g_press(is:ie,js:je,ks:ke+1)=0.0
183    do k=ks,ke
184       do j=js,je
185          do i=is,ie
186             g_press(i,j,k+1)=g_press(i,j,k+1)+0.5*grid%xa%p(i,j,k)
187             g_press(i,j,k )=g_press(i,j,k )+0.5*grid%xa%p(i,j,k)
188             grid%g_moist(i,j,k,P_G_QV)=grid%g_moist(i,j,k,P_G_QV)-(grid%em_mu_2(i,j)+grid%em_mub(i,j))*g_press(i,j,k)*grid%em_dn(k)
189             grid%g_mu_2(i,j)=grid%g_mu_2(i,j)-g_press(i,j,k)*(1.0+grid%moist(i,j,k,P_QV))*grid%em_dn(k)
190             g_press(i,j,k+1)=g_press(i,j,k+1)+g_press(i,j,k)
191          end do
192       end do
193    end do
194 
195    grid%xa%p(is:ie, js:je, ks:ke) = 0.0
196 
197    !---------------------------------------------------------------------------
198    ! [2.0] COMPUTE psfc increments from mu-increments
199    !---------------------------------------------------------------------------
200 
201    do j=js,je
202       do i=is,ie
203          sdmd=0.0
204          s1md=0.0
205          do k=ks,ke
206             s1md=s1md+(1.0+grid%moist(i,j,k,P_QV))*grid%em_dnw(k)
207          end do
208          grid%g_mu_2(i,j)=grid%g_mu_2(i,j)-grid%xa%psfc(i,j)*s1md
209          sdmd=sdmd-grid%xb%psac(i,j)*grid%xa%psfc(i,j)
210          do k=ks,ke
211             grid%g_moist(i,j,k,P_G_QV)=grid%g_moist(i,j,k,P_G_QV)+sdmd*grid%em_dnw(k)
212          end do
213       end do
214    end do
215 
216    grid%xa%psfc(is:ie, js:je) = 0.0
217 
218    !---------------------------------------------------------------------------
219    ! [1.0] Get the specific humidity increments from mixing ratio increments
220    !---------------------------------------------------------------------------
221    do k=ks,ke
222       do j=js,je
223          do i=is,ie
224             grid%g_moist(i,j,k,P_G_QV)=grid%g_moist(i,j,k,P_G_QV)+grid%xa%q(i,j,k)* &
225                (1.0-grid%xb%q(i,j,k))**2
226          end do
227       end do
228    end do
229 
230    grid%xa%q(is:ie, js:je, ks:ke) = 0.0
231 
232    ndynopt      = grid%dyn_opt
233    grid%dyn_opt = DYN_EM_TL
234    call nl_set_dyn_opt (1 , DYN_EM_TL)
235 
236    call da_med_initialdata_output(grid , config_flags, filnam)
237 
238    grid%dyn_opt = ndynopt
239    call nl_set_dyn_opt (1 , DYN_EM)
240 
241    if (trace_use) call da_trace_exit("da_transfer_wrftltoxa_adj")
242 
243 end subroutine da_transfer_wrftltoxa_adj
244 
245