wrf svn trunk commit r4103
[wrffire.git] / wrfv2_fire / var / da / da_vtox_transforms / da_transform_vptox_adj.inc
blobcb83d2b3c867f58c64ccf105780c1af0290892cf
1 subroutine da_transform_vptox_adj(grid, vp, be, ep)
3    !--------------------------------------------------------------------------
4    ! Purpose: Adjoint for Physical transform of 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    type (vp_type), intent(inout)        :: vp  ! CV on grid structure.
13    type (ep_type), intent(in)           :: ep  ! Ensemble perturbation.
14    type (be_type), intent(in), optional :: be  ! Background errors.
16    integer :: i, k, j, ij, k1              ! Loop counters.
18    if (trace_use) call da_trace_entry("da_transform_vptox_adj")
20    !---------------------------------------------------------------------------
21    !  [4] Add flow-dependent increments in model space (grid%xa):
22    !---------------------------------------------------------------------------
23       
24    if (be % ne > 0 .and. alphacv_method == alphacv_method_xa) then
25       call da_add_flow_dependence_xa_adj(be % ne, ep, grid%xa, vp)
26    end if
28    !--------------------------------------------------------------------------
29    ! [3] Transform to model variable space:
30    !--------------------------------------------------------------------------
32    if ((use_radarobs .and. use_radar_rf) .or. (use_rad .and. crtm_cloud) ) then
33       ! Pseudo RH --> Total water mixing ratio:
34       vp%v4(its:ite,jts:jte,kts:kte)  = vp%v4(its:ite,jts:jte,kts:kte) &
35          + grid%xa%qt(its:ite,jts:jte,kts:kte) * grid%xb%qs(its:ite,jts:jte,kts:kte)
36    else 
37       ! Pseudo RH --> Water vapor mixing ratio:
38       !$OMP PARALLEL DO &
39       !$OMP PRIVATE ( ij )
40       do ij = 1 , grid%num_tiles
41          do k = kts,kte
42             do j = grid%j_start(ij),grid%j_end(ij)
43                do i =  its, ite
44                   vp%v4(i,j,k)  = vp%v4(i,j,k) + grid%xa%q(i,j,k) * grid%xb%qs(i,j,k)   
45                end do
46             end do
47          end do
48       end do
49       !$OMP END PARALLEL DO
50    end if
52 #ifdef A2C
53   if ((fg_format==fg_format_wrf_arw_regional  .or. &
54        fg_format==fg_format_wrf_arw_global  ) .and. ide == ipe ) then
55      ipe = ipe + 1
56      ide = ide + 1
57   end if
59   if ((fg_format==fg_format_wrf_arw_regional  .or. &
60        fg_format==fg_format_wrf_arw_global  ) .and. jde == jpe ) then
61      jpe = jpe + 1
62      jde = jde + 1
63   end if
64 #endif
65 #ifdef DM_PARALLEL
66 #include "HALO_PSICHI_UV_ADJ.inc"
67 #endif
69 #ifdef A2C
70   if ((fg_format==fg_format_wrf_arw_regional  .or. &
71        fg_format==fg_format_wrf_arw_global  ) .and. ide == ipe ) then
72      ipe = ipe - 1
73      ide = ide - 1
74   end if
76   if ((fg_format==fg_format_wrf_arw_regional  .or. &
77        fg_format==fg_format_wrf_arw_global  ) .and. jde == jpe ) then
78      jpe = jpe - 1
79      jde = jde - 1
80   end if
81 #endif
84    ! Transform psi and chi to u and v:
86    call da_psichi_to_uv_adj(grid%xa % u, grid%xa % v, grid%xb % coefx, grid%xb % coefy, vp % v1, vp % v2)
88    !--------------------------------------------------------------------------
89    ! [2] Impose statistical balance constraints:
90    !--------------------------------------------------------------------------
92    !$OMP PARALLEL DO &
93    !$OMP PRIVATE ( ij, k, j, k1, i )
94    do ij = 1 , grid%num_tiles
96    ! Surface Pressure
97    do k= kts,kte
98       do j= grid%j_start(ij),grid%j_end(ij)
99          do i =  its, ite
100             vp%v1(i,j,k) = vp%v1(i,j,k) + be%reg_ps(j,k)*grid%xa%psfc(i,j)
101          end do
102       end do
103    end do
104    do j= grid%j_start(ij),grid%j_end(ij)
105       do i =  its, ite
106          vp%v5(i,j,1) = grid%xa%psfc(i,j) 
107       end do
108    end do
110    ! Temperature
111    do k1 = kts,kte
112       do k = kts,kte
113          do j = grid%j_start(ij),grid%j_end(ij)
114             do i =  its, ite
115                vp%v1(i,j,k1) = vp%v1(i,j,k1) + be%reg_t(j,k,k1)*grid%xa%t(i,j,k)
116             end do
117          end do
118       end do
119    end do
120    do k = kts,kte
121       do j = grid%j_start(ij),grid%j_end(ij)
122          do i =  its, ite
123             vp%v3(i,j,k) = grid%xa%t(i,j,k)
124          end do
125       end do
126    end do
128    ! Chi
129    do k = kts,kte
130       do j = grid%j_start(ij),grid%j_end(ij)
131          do i =  its, ite
132             vp%v1(i,j,k) = vp%v1(i,j,k) + be%reg_chi(j,k)*vp%v2(i,j,k)
133          end do
134        end do
135    end do
137    enddo
138    !$OMP END PARALLEL DO
140    !---------------------------------------------------------------------------
141    !  [1] Add flow-dependent increments in control variable space (vp):
142    !---------------------------------------------------------------------------
143    
144    if (be % ne > 0 .and. alphacv_method == alphacv_method_vp) then
145       call da_add_flow_dependence_vp_adj(be % ne, ep, vp)
146    end if
148    if (trace_use) call da_trace_exit("da_transform_vptox_adj")
150 end subroutine da_transform_vptox_adj