wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_vtox_transforms / da_transform_vptox.inc
blob38722e028034d9db6369d38f6670119a3f2c0bd3
1 subroutine da_transform_vptox(grid, vp, be, ep)
3    !-----------------------------------------------------------------------
4    ! Purpose: Physical transform of analysis increment variables.
5    !    Updated for Analysis on Arakawa-C grid
6    !    Author: Syed RH Rizvi,  MMM/ESSL/NCAR,  Date: 10/22/2008
7    !-----------------------------------------------------------------------
9    implicit none
11    type (domain), intent(inout)         :: grid
12    
13    type (vp_type), intent(inout)        :: vp  ! CV on grid structure.
14    type (be_type), intent(in), optional :: be  ! Background errors.
15    type (ep_type), intent(in), optional :: ep  ! Ensemble perturbations.
17    integer :: i, k, j, k1, ij            ! Loop counters.
19    if (trace_use) call da_trace_entry("da_transform_vptox") 
21    !---------------------------------------------------------------------------
22    !  [1] Add flow-dependent increments in control variable space (vp):
23    !---------------------------------------------------------------------------
25    if (be % ne > 0 .and. alphacv_method == alphacv_method_vp) then
26       call da_add_flow_dependence_vp(be % ne, ep, vp, its,ite, jts,jte, kts,kte)
27    end if
29    !--------------------------------------------------------------------------
30    ! [2] Impose statistical balance constraints:
31    !--------------------------------------------------------------------------
33    !$OMP PARALLEL DO &
34    !$OMP PRIVATE ( ij, k1, k, j, i)
35    do ij = 1 , grid%num_tiles
37    ! Chi:
38    do k = kts, kte
39       do j = grid%j_start(ij), grid%j_end(ij)
40          do i = its, ite
41             vp%v2(i,j,k) = vp%v2(i,j,k) + be%reg_chi(j,k)* vp%v1(i,j,k)
42          end do
43       end do
44    end do
45   
46    ! Temperature:
47    do k = kts, kte
48       do j = grid%j_start(ij), grid%j_end(ij)
49          do i = its, ite
50             grid%xa%t(i,j,k) = vp%v3(i,j,k)
51          end do
52       end do
53    end do
56    do k1 = kts, kte
57       do k = kts, kte
58          do j = grid%j_start(ij), grid%j_end(ij)
59             do i = its, ite
60                grid%xa%t(i,j,k) = grid%xa%t(i,j,k) + be%reg_t(j,k,k1)*vp%v1(i,j,k1)
61             end do
62          end do
63       end do
64    end do
66    ! Surface Pressure
67    do j = grid%j_start(ij), grid%j_end(ij)
68       do i = its, ite
69          grid%xa%psfc(i,j) = vp%v5(i,j,1) 
70       end do
71    end do
73    do k = kts,kte
74       do j = grid%j_start(ij), grid%j_end(ij)
75          do i = its, ite
76             grid%xa%psfc(i,j) = grid%xa%psfc(i,j) + be%reg_ps(j,k)*vp%v1(i,j,k)
77          end do
78       end do
79    end do
81    end do
82    !$OMP END PARALLEL DO
83    !--------------------------------------------------------------------------
84    ! [3] Transform to model variable space:
85    !--------------------------------------------------------------------------
86   
87 #ifdef A2C
88   if ((fg_format==fg_format_wrf_arw_regional  .or. &
89        fg_format==fg_format_wrf_arw_global  ) .and. ide == ipe ) then
90      ipe = ipe + 1
91      ide = ide + 1
92   end if
94   if ((fg_format==fg_format_wrf_arw_regional  .or. &
95        fg_format==fg_format_wrf_arw_global  ) .and. jde == jpe ) then
96      jpe = jpe + 1
97      jde = jde + 1
98   end if
99 #endif
101 #ifdef DM_PARALLEL
102 #include "HALO_PSICHI_UV.inc"
103 #endif
105 #ifdef A2C
106   if ((fg_format==fg_format_wrf_arw_regional  .or. &
107        fg_format==fg_format_wrf_arw_global  ) .and. ide == ipe ) then
108      ipe = ipe - 1
109      ide = ide - 1
110   end if
112   if ((fg_format==fg_format_wrf_arw_regional  .or. &
113        fg_format==fg_format_wrf_arw_global  ) .and. jde == jpe ) then
114      jpe = jpe - 1
115      jde = jde - 1
116   end if
117 #endif
119    ! Psi and chi to u and v:
121    call da_psichi_to_uv(vp % v1, vp % v2, grid%xb % coefx, &
122       grid%xb % coefy , grid%xa % u, grid%xa % v)
124    if ((use_radarobs .and. use_radar_rf) .or. (use_rad .and. crtm_cloud)) then
125       ! Pseudo RH --> Total water mixing ratio:
126       grid%xa % qt(its:ite,jts:jte,kts:kte) = vp%v4(its:ite,jts:jte,kts:kte) * &
127          grid%xb%qs(its:ite,jts:jte,kts:kte)
128    else
129       ! Pseudo RH --> Water vapor mixing ratio:
130       !$OMP PARALLEL DO &
131       !$OMP PRIVATE ( ij, i, j, k )
132       do ij = 1 , grid%num_tiles
133          do k = kts, kte
134             do j = grid%j_start(ij), grid%j_end(ij)
135                do i = its, ite
136                   grid%xa % q(i,j,k) =  vp%v4(i,j,k) * grid%xb%qs(i,j,k)
137                enddo
138             enddo
139          enddo
140       enddo
141       !$OMP END PARALLEL DO
142    end if
144    !---------------------------------------------------------------------------
145    !  [4] Add flow-dependent increments in model space (grid%xa):
146    !---------------------------------------------------------------------------
148    if (be % ne > 0 .and. alphacv_method == alphacv_method_xa) then
149       call da_add_flow_dependence_xa(grid, be % ne, ep, vp)
150    end if
152    if (trace_use) call da_trace_exit("da_transform_vptox") 
154 end subroutine da_transform_vptox