da_transform_xtogpsref_adj.inc
References to this file elsewhere.
1 subroutine da_transform_xtogpsref_adj(xa, xb, xp)
2
3 !-----------------------------------------------------------------------
4 ! Purpose: TBD
5 ! input: xb%q, xb%p, xb%t, xa%ref, and xp
6 ! output: xa%q, xa%p, xa%t
7 !-----------------------------------------------------------------------
8
9 implicit none
10
11 type (x_type), intent(inout) :: xa ! gridded analysis increment.
12 type (xb_type), intent(in) :: xb ! first guess state.
13 type (xpose_type), intent(in) :: xp ! domain decomposition vars.
14
15 integer :: i, j, k, is, ie, js, je, ks, ke
16 real :: partone, parttwo, dividnd
17 real :: partone9,parttwo9,dividnd9
18
19 if (trace_use) call da_trace_entry("da_transform_xtogpsref_adj")
20
21 ! 1.0 Get the range for a tile
22
23 is = xp%its; ie = xp%ite
24 js = xp%jts; je = xp%jte
25 ks = xp%kts; ke = xp%kte
26
27 if (Testing_WRFVAR) then
28 is = xb%its-1
29 js = xb%jts-1
30
31 ie = xb%ite+1
32 je = xb%jte+1
33
34 if (is < xb%ids) is = xb%ids
35 if (js < xb%jds) js = xb%jds
36
37 if (ie > xb%ide) ie = xb%ide
38 if (je > xb%jde) je = xb%jde
39 end if
40
41 !------------------------------------------------------------------------------
42 ! [2.0] Calculate the adjoint for GPS refractivity:
43 !------------------------------------------------------------------------------
44
45 do k = ks, ke
46 do j = js, je
47 do i = is, ie
48
49 ! Note: p in Pascal, q is the specific humidity
50
51 ! 2.1 basic state
52
53 partone9 = 0.776*xb%p(i,j,k)/xb%t(i,j,k)
54 dividnd9 = xb%t(i,j,k)*(0.622+0.378*xb%q(i,j,k))
55 parttwo9 = 1.0+coeff*xb%q(i,j,k)/dividnd9
56
57 ! 2.2 Adjoint of partone and parttwo
58
59 partone = xa%ref(i,j,k) * parttwo9
60 parttwo = xa%ref(i,j,k) * partone9
61
62 ! 2.3 Adjoint of q and dividnd
63
64 xa%q(i,j,k) = xa%q(i,j,k) + coeff*parttwo/dividnd9
65 dividnd=-coeff*xb%q(i,j,k)*parttwo/(dividnd9*dividnd9)
66
67 ! 2.4 Adjoint of t, q, and p
68
69 xa%t(i,j,k) = xa%t(i,j,k) + dividnd*(0.622+0.378*xb%q(i,j,k))
70 xa%q(i,j,k) = xa%q(i,j,k) + xb%t(i,j,k)*0.378*dividnd
71
72 xa%p(i,j,k) = xa%p(i,j,k) + 0.776*partone/xb%t(i,j,k)
73 xa%t(i,j,k) = xa%t(i,j,k) - &
74 0.776*xb%p(i,j,k)*partone/(xb%t(i,j,k)*xb%t(i,j,k))
75
76 ! [3.0] Clean the xa%ref
77
78 xa%ref(i,j,k) = 0.0
79
80 end do
81 end do
82 end do
83
84 if (trace_use) call da_trace_exit("da_transform_xtogpsref_adj")
85
86 end subroutine da_transform_xtogpsref_adj
87
88