da_transform_vptox_adj.inc

References to this file elsewhere.
1 subroutine da_transform_vptox_adj(grid, vp, be, ep)
2 
3    !--------------------------------------------------------------------------
4    ! Purpose: Adjoint for Physical transform of variables 
5    !--------------------------------------------------------------------------
6 
7    implicit none
8 
9    type (domain),  intent(inout)        :: grid
10    type (vp_type), intent(inout)        :: vp  ! CV on grid structure.
11    type (ep_type), intent(in)           :: ep  ! Ensemble perturbation.
12    type (be_type), intent(in), optional :: be  ! Background errors.
13 
14    integer :: k, j, k1              ! Loop counters.
15 
16    if (trace_use) call da_trace_entry("da_transform_vptox_adj")
17 
18    !---------------------------------------------------------------------------
19    !  [4] Add flow-dependent increments in model space (grid%xa):
20    !---------------------------------------------------------------------------
21       
22    if (be % ne > 0 .and. alphacv_method == 2) then
23       call da_add_flow_dependence_xa_adj(be % ne, ep, grid%xa, vp)
24    end if
25 
26    !--------------------------------------------------------------------------
27    ! [3] Transform to model variable space:
28    !--------------------------------------------------------------------------
29 
30    if (use_radarobs .and. use_radar_rf) then
31       ! Pseudo RH --> Total water mixing ratio:
32       vp%v4(its:ite,jts:jte,kts:kte)  = vp%v4(its:ite,jts:jte,kts:kte) &
33          + grid%xa%qt(its:ite,jts:jte,kts:kte) * grid%xb%qs(its:ite,jts:jte,kts:kte)
34    else 
35       ! Pseudo RH --> Water vapor mixing ratio:
36       vp%v4(its:ite,jts:jte,kts:kte)  = vp%v4(its:ite,jts:jte,kts:kte) &
37          + grid%xa%q (its:ite,jts:jte,kts:kte) * grid%xb%qs(its:ite,jts:jte,kts:kte)
38    end if
39 
40 #ifdef DM_PARALLEL
41 #include "HALO_PSICHI_UV_ADJ.inc"
42 #endif
43 
44    ! Transform psi and chi to u and v:
45    call da_psichi_to_uv_adj(grid%xa % u, grid%xa % v, grid%xb % coefx, grid%xb % coefy, vp % v1, vp % v2)
46 
47    !--------------------------------------------------------------------------
48    ! [2] Impose statistical balance constraints:
49    !--------------------------------------------------------------------------
50 
51    ! Surface Pressure
52    do k=kts,kte
53       do j=jts,jte
54          vp%v1(its:ite,j,k) = vp%v1(its:ite,j,k) + be%reg_ps(j,k)*grid%xa%psfc(its:ite,j)
55       end do
56    end do
57    vp%v5(its:ite,jts:jte,1) = grid%xa%psfc(its:ite,jts:jte) 
58 
59    ! Temperature
60    do k1 = kts,kte
61       do k = kts,kte
62          do j = jts,jte
63             vp%v1(its:ite,j,k1) = vp%v1(its:ite,j,k1) + be%reg_t(j,k,k1)*grid%xa%t(its:ite,j,k)
64          end do
65       end do
66    end do
67    vp%v3(its:ite,jts:jte,kts:kte) = grid%xa%t(its:ite,jts:jte,kts:kte)
68 
69    ! Chi
70    do k = kts,kte
71       do j = jts,jte
72          vp%v1(its:ite,j,k) = vp%v1(its:ite,j,k) + be%reg_chi(j,k)*vp%v2(its:ite,j,k)
73        end do
74    end do
75 
76    !---------------------------------------------------------------------------
77    !  [1] Add flow-dependent increments in control variable space (vp):
78    !---------------------------------------------------------------------------
79    
80    if (be % ne > 0 .and. alphacv_method == alphacv_method_vp) then
81       call da_add_flow_dependence_vp_adj(be % ne, ep, vp)
82    end if
83 
84    if (trace_use) call da_trace_exit("da_transform_vptox_adj")
85 
86 end subroutine da_transform_vptox_adj
87 
88