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