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