da_moist_phys_lin.inc

References to this file elsewhere.
1 subroutine da_moist_phys_lin(xb, xa, xp)
2 
3    !---------------------------------------------------------------------------
4    !  Purpose: Partition of the hydrometeors via the moist explicit scheme.
5    !           A warm rain process is used in this subroutine.
6    !           This is the tangent linear code of the scheme.
7    !
8    !  Method: The warm rain process is according to Hsie and Anthes (1984)
9    !          and Dudhia (1989)
10    !
11    !  Assumptions: 1) Model level stored top down.
12    !---------------------------------------------------------------------------
13 
14    implicit none
15 
16    type (xb_type), intent(in)    :: xb           ! First guess structure.
17    type (x_type), intent(inout)  :: xa           ! increment structure.
18    type (xpose_type), intent(inout) :: xp        ! Dimensions and xpose buffers.
19 
20    real, dimension(ims:ime,jms:jme,kms:kme) :: T_OLD,T_NEW
21    real, dimension(ims:ime,jms:jme,kms:kme) :: Q_OLD,Q_NEW
22    real, dimension(ims:ime,jms:jme,kms:kme) :: QCW_OLD,QCW_NEW
23    real, dimension(ims:ime,jms:jme,kms:kme) :: QRN_OLD,QRN_NEW
24 
25    real, dimension(kms:kme)                 :: EES, QVSS
26    real, dimension(kms:kme)                 :: EES9, QVSS9
27 
28    real, dimension(kms:kme)                   :: DT
29    real, dimension(kms:kme)                   :: QVT,QCT,QRT,TTT
30    real, dimension(kms:kme)                   :: QVT9,QCT9,QRT9,TTT9
31    real, dimension(kms:kme) :: SCR2,SCR3,SCR4,SCR5,SCR6
32    real, dimension(kms:kme) :: DUM31
33    real, dimension(kms:kme) :: PRA,PRC,PRD,PRE
34    real, dimension(kms:kme) :: SCR31,SCR42,SCR71
35    real, dimension(kms:kme) :: DUM112,DUM113,DUM211,DUM411
36 
37    real, dimension(kms:kme) :: SCR29,SCR39,SCR49,SCR59,SCR69
38    real, dimension(kms:kme) :: DUM319
39    real, dimension(kms:kme) :: PRA9,PRC9,PRD9,PRE9
40    real, dimension(kms:kme) :: SCR319,SCR429,SCR719
41    real, dimension(kms:kme) :: DUM1129,DUM1139,DUM2119,DUM4119
42    real, dimension(kms:kme) :: TMP
43 
44    integer :: i, j, k
45 
46    do k=kts,kte
47       do j=jts,jte
48          do i=its,ite
49             T_OLD(i,j,k) = xa % t (i,j,k)
50             Q_OLD(i,j,k) = xa % q (i,j,k)
51             QCW_OLD(i,j,k) = xa % qcw (i,j,k)
52             QRN_OLD(i,j,k) = xa % qrn (i,j,k)
53          end do
54       end do
55    end do
56 
57    do j=jts,jte
58       do i=its,ite
59 
60       !  Preparation
61 
62          do K=kts,kte
63             xa % q(i,j,k) =xa % qt(i,j,k) - xa % qcw(i,j,k) - xa %qrn(i,j,k)
64          end do
65 
66          do K=kts,kte
67             DT(k) = xb%delt(i,j,k)
68          end do
69 
70          do K=kts,kte
71 
72             if (DT(k) <= 0.) cycle
73 
74             if ( xb%t(I,J,K) > TO )then
75                EES(K)=SVP1*EXP(SVP2*(xb%t(I,J,K)-SVPT0)/(xb%t(I,J,K)-SVP3))
76                EES9(K)=EES(K)*SVP2*(SVPT0-SVP3)/((xb%t(I,J,K)-SVP3)*  &
77                        (xb%t(I,J,K)-SVP3))*xa%t(I,J,K)
78             else
79                EES(K)=.611*EXP(22.514-6.15E3/xb%t(I,J,K))
80                EES9(K)=EES(K)*6.15E3/(xb%t(I,J,K)*xb%t(I,J,K))*xa%t(I,J,K)
81             end if
82 
83             TMP(K)=622./((xb%p(I,J,K)-EES(K))**2)
84             QVSS9(K)=TMP(K)*xb%p(I,J,K)*EES9(K)  &
85                   -TMP(K)*EES(K)*xa%p(I,J,K)
86             QVSS(K)=622.*EES(K)/(xb%p(I,J,K)-EES(K))
87 
88 
89             SCR49(K)=xa%q(I,J,K)/QVSS(K)-xb%q(I,J,K)/QVSS(K)**2*QVSS9(K)
90             SCR4(K)=xb%q(I,J,K)/QVSS(K)
91          end do
92 
93          do K=kts,kte
94 
95             if (DT(k) <= 0.) cycle
96 
97             if (xb%qcw(I,J,K) > 0.) then
98                SCR29(K)=xa%qcw(I,J,K)
99                SCR2(K)=xb%qcw(I,J,K)
100             else
101                SCR29(K)=0.
102                SCR2(K)=0.
103             end if
104             if (xb%qrn(I,J,K) > 1.E-25) then
105                SCR39(K)=xa%qrn(I,J,K)
106                SCR3(K)=xb%qrn(I,J,K)
107             else
108                SCR39(K)=0.
109                SCR3(K)=1.E-25
110             end if
111             SCR59(K)=xa%q(I,J,K)/SCR4(K)-xb%q(I,J,K)/SCR4(K)**2*SCR49(K)
112             SCR5(K)=xb%q(I,J,K)/SCR4(K)
113 
114             SCR69(K)=xa%p(I,J,K)/(gas_constant*xb%t(I,J,K))-xb%p(I,J,K)/  &
115                      (gas_constant*xb%t(I,J,K)**2)*xa%t(I,J,K)
116             SCR6(K)=xb%p(I,J,K)/(gas_constant*xb%t(I,J,K))
117 
118             DUM319(K)=-XLV1*xa%t(I,J,K) 
119             DUM31(K)=3.1484E6-XLV1*xb%t(I,J,K)
120 
121          end do
122 
123          call da_autoc_lin(DT,SCR2,PRC,SCR29,PRC9,kts,kte,kms,kme)
124 
125          call da_accre_lin(DT,SCR2,SCR3,PRA,SCR29,SCR39,PRA9,kts,kte,kms,kme)
126 
127          call da_evapo_lin(DT,SCR3,SCR5,xb%q(I,J,:),PRE,SCR6,  &
128                            SCR39,SCR59,xa%q(I,J,:),PRE9,SCR69, &
129                            kts,kte,kms,kme)
130 
131          do K=kts, kte
132 
133             if (DT(k) <= 0.) cycle
134 
135             !  Readjust
136 
137             DUM1129(K)=(PRC9(k)+PRA9(k))*DT(k)
138             DUM112(K)=(PRC(k)+PRA(k))*DT(k)
139             if (DUM112(K) > SCR2(k)) then
140                PRC9(K)=SCR29(K)*PRC(K)/DUM112(K)  &
141                       +PRC9(K)*SCR2(K)/DUM112(K)  &
142                       -SCR2(K)*PRC(K)/(DUM112(K)*DUM112(K))*DUM1129(K)
143                PRC(k)=SCR2(K)*PRC(K)/DUM112(K)
144                PRA9(K)=SCR29(K)*PRA(K)/DUM112(K)  &
145                       +PRA9(K)*SCR2(K)/DUM112(K)  &
146                       -SCR2(K)*PRA(K)/(DUM112(K)*DUM112(K))*DUM1129(K)
147                PRA(k)=SCR2(K)*PRA(K)/DUM112(K)
148             end if
149             QVT9(K)=-PRE9(K)
150             QVT(K)=-PRE(K)
151             QCT9(K)=-PRC9(K)-PRA9(K)
152             QCT(K)=-PRC(K)-PRA(K)
153             QRT9(K)=PRC9(K)+PRA9(K)+PRE9(K)
154             QRT(K)=PRC(K)+PRA(K)+PRE(K)
155             if (xb%t(I,J,K).GT.TO)then
156                DUM4119(K)=DUM319(K)
157                DUM411(K)=DUM31(K)
158             else
159                DUM4119(K)=0.
160                DUM411(K)=XLS
161             end if
162             PRD9(K)=cp*0.887*xa%q(I,J,K)
163             PRD(K)=cp*(1.+0.887*xb%q(I,J,K))
164             TTT9(K)=-DUM4119(K)*QVT(K)/PRD(K)  &
165                    -QVT9(K)*DUM411(K)/PRD(K)  &
166                    +DUM411(K)*QVT(K)/(PRD(K)*PRD(K))*PRD9(K)
167             TTT(K)=-DUM411(K)*QVT(K)/PRD(K)
168 
169             DUM1139(K)=xa%q(I,J,K)+DT(K)*QVT9(K)
170             DUM113(K)=xb%q(I,J,K)+DT(K)*QVT(K)
171             if (DUM113(K) > 1.E-12 ) then
172                SCR429(K)=DUM1139(K)
173                SCR42(K)=DUM113(K)
174             else
175                SCR429(K)=0.
176                SCR42(K)=1.E-12
177             end if
178             DUM2119(K)=xa%qcw(I,J,K)+QCT9(K)*DT(k)
179             DUM211(K)=xb%qcw(I,J,K)+QCT(K)*DT(k)
180             if (DUM211(K) > 0.) then
181                SCR319(K)=DUM2119(K)
182                SCR31(K)=DUM211(K)
183             else
184                SCR319(K)=0.
185                SCR31(K)=0.
186             end if
187             SCR719(K)=xa%t(I,J,K)+TTT9(K)*DT(k)
188             SCR71(K)=xb%t(I,J,K)+TTT(K)*DT(k)
189          end do
190 
191          call da_condens_lin(DT,SCR31,SCR42,SCR71,DUM31,PRD,         &
192                              QVT,QCT,QRT,TTT,                        &
193                              xb%p(I,J,:),xb%t(I,J,:),xb%q(I,J,:),    &
194                              xb%qcw(I,J,:),xb%qrn(I,J,:),            &
195                              SCR319,SCR429,SCR719,DUM319,PRD9,       &
196                              QVT9,QCT9,QRT9,TTT9,                    &
197                              xa%p(I,J,:),xa%t(I,J,:),xa%q(I,J,:),    &
198                              xa%qcw(I,J,:),xa%qrn(I,J,:),kts,kte)
199       end do
200    end do
201 
202    do k=kds,kde
203       do j=jts,jte
204          do i=its,ite
205             T_NEW(i,j,k) = xa % t (i,j,k) - T_OLD(i,j,k)
206             Q_NEW(i,j,k) = xa % q (i,j,k) - Q_OLD(i,j,k)
207             QCW_NEW(i,j,k) = xa % qcw (i,j,k) - QCW_OLD(i,j,k)
208             QRN_NEW(i,j,k) = xa % qrn (i,j,k) - QRN_OLD(i,j,k)
209          end do
210       end do
211    end do
212 
213    call da_filter(t_new, xp)
214    call da_filter(q_new, xp)
215    call da_filter(qcw_new, xp)
216    call da_filter(qrn_new, xp)
217 
218    do k=kds,kde
219       do j=jts,jte
220          do i=its,ite
221             xa % t (i,j,k) = T_NEW(i,j,k) + T_OLD(i,j,k)
222             xa % q (i,j,k) = Q_NEW(i,j,k) + Q_OLD(i,j,k)
223             xa % qcw (i,j,k) = QCW_NEW(i,j,k) + QCW_OLD(i,j,k)
224             xa % qrn (i,j,k) = QRN_NEW(i,j,k) + QRN_OLD(i,j,k)
225          end do
226       end do
227    end do
228 
229 end subroutine da_moist_phys_lin
230 
231