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