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