1 |
module cv30_unsat_m |
2 |
|
3 |
implicit none |
4 |
|
5 |
contains |
6 |
|
7 |
SUBROUTINE cv30_unsat(icb, inb, t, q, qs, gz, u, v, p, ph, th, tv, lv, cpn, & |
8 |
ep, sigp, clw, m, ment, elij, delt, plcl, mp, qp, up, vp, wt, water, & |
9 |
evap, b) |
10 |
|
11 |
use cv30_param_m, only: nl, sigd |
12 |
use cv_thermo_m, only: cpd, ginv, grav |
13 |
USE dimphy, ONLY: klon, klev |
14 |
|
15 |
integer, intent(in):: icb(:) ! (ncum) |
16 |
|
17 |
integer, intent(in):: inb(:) ! (ncum) |
18 |
! first model level above the level of neutral buoyancy of the |
19 |
! parcel (1 <= inb <= nl - 1) |
20 |
|
21 |
real, intent(in):: t(:, :), q(:, :), qs(:, :) ! (klon, klev) |
22 |
real, intent(in):: gz(:, :) ! (klon, klev) |
23 |
real, intent(in):: u(:, :), v(:, :) ! (klon, klev) |
24 |
real, intent(in):: p(klon, klev), ph(klon, klev + 1) |
25 |
real, intent(in):: th(klon, klev) |
26 |
real, intent(in):: tv(klon, klev) |
27 |
real, intent(in):: lv(klon, klev) |
28 |
real, intent(in):: cpn(klon, klev) |
29 |
real, intent(in):: ep(:, :), sigp(:, :), clw(:, :) ! (ncum, klev) |
30 |
real, intent(in):: m(:, :) ! (ncum, klev) |
31 |
real, intent(in):: ment(:, :, :) ! (ncum, klev, klev) |
32 |
real, intent(in):: elij(:, :, :) ! (ncum, klev, klev) |
33 |
real, intent(in):: delt |
34 |
real, intent(in):: plcl(klon) |
35 |
|
36 |
! outputs: |
37 |
real, intent(out):: mp(klon, klev) |
38 |
real, intent(out):: qp(:, :), up(:, :), vp(:, :) ! (ncum, nl) |
39 |
real wt(klon, klev), water(klon, klev), evap(klon, klev) |
40 |
real, intent(out):: b(:, :) ! (ncum, nl - 1) |
41 |
|
42 |
! Local: |
43 |
integer ncum |
44 |
integer i, j, il, num1 |
45 |
real tinv, delti |
46 |
real awat, afac, afac1, afac2, bfac |
47 |
real pr1, pr2, sigt, b6, c6, revap, tevap, delth |
48 |
real amfac, amp2, xf, tf, fac2, ur, sru, fac, d, af, bf |
49 |
real ampmax |
50 |
real lvcp(klon, klev) |
51 |
real wdtrain(size(icb)) |
52 |
logical lwork(size(icb)) |
53 |
|
54 |
!------------------------------------------------------ |
55 |
|
56 |
ncum = size(icb) |
57 |
delti = 1. / delt |
58 |
tinv = 1. / 3. |
59 |
mp = 0. |
60 |
b = 0. |
61 |
|
62 |
do i = 1, nl |
63 |
do il = 1, ncum |
64 |
qp(il, i) = q(il, i) |
65 |
up(il, i) = u(il, i) |
66 |
vp(il, i) = v(il, i) |
67 |
wt(il, i) = 0.001 |
68 |
water(il, i) = 0. |
69 |
evap(il, i) = 0. |
70 |
lvcp(il, i) = lv(il, i) / cpn(il, i) |
71 |
enddo |
72 |
enddo |
73 |
|
74 |
! check whether ep(inb) = 0, if so, skip precipitating |
75 |
! downdraft calculation |
76 |
forall (il = 1:ncum) lwork(il) = ep(il, inb(il)) >= 1e-4 |
77 |
|
78 |
wdtrain = 0. |
79 |
|
80 |
downdraft_loop: DO i = nl - 1, 1, - 1 |
81 |
num1 = 0 |
82 |
|
83 |
do il = 1, ncum |
84 |
if (i <= inb(il) .and. lwork(il)) num1 = num1 + 1 |
85 |
enddo |
86 |
|
87 |
if (num1 > 0) then |
88 |
! integrate liquid water equation to find condensed water |
89 |
! and condensed water flux |
90 |
|
91 |
! calculate detrained precipitation |
92 |
|
93 |
do il = 1, ncum |
94 |
if (i <= inb(il) .and. lwork(il)) then |
95 |
wdtrain(il) = grav * ep(il, i) * m(il, i) * clw(il, i) |
96 |
endif |
97 |
enddo |
98 |
|
99 |
if (i > 1) then |
100 |
do j = 1, i - 1 |
101 |
do il = 1, ncum |
102 |
if (i <= inb(il) .and. lwork(il)) then |
103 |
awat = elij(il, j, i) - (1. - ep(il, i)) * clw(il, i) |
104 |
awat = max(awat, 0.) |
105 |
wdtrain(il) = wdtrain(il) + grav * awat * ment(il, j, i) |
106 |
endif |
107 |
enddo |
108 |
end do |
109 |
endif |
110 |
|
111 |
! find rain water and evaporation using provisional |
112 |
! estimates of qp(i) and qp(i - 1) |
113 |
|
114 |
do il = 1, ncum |
115 |
if (i <= inb(il) .and. lwork(il)) then |
116 |
wt(il, i) = 45. |
117 |
|
118 |
if (i < inb(il)) then |
119 |
qp(il, i) = qp(il, i + 1) + (cpd * (t(il, i + 1) & |
120 |
- t(il, i)) + gz(il, i + 1) - gz(il, i)) / lv(il, i) |
121 |
qp(il, i) = 0.5 * (qp(il, i) + q(il, i)) |
122 |
endif |
123 |
|
124 |
qp(il, i) = max(qp(il, i), 0.) |
125 |
qp(il, i) = min(qp(il, i), qs(il, i)) |
126 |
qp(il, inb(il)) = q(il, inb(il)) |
127 |
|
128 |
if (i == 1) then |
129 |
afac = p(il, 1) * (qs(il, 1) - qp(il, 1)) & |
130 |
/ (1e4 + 2000. * p(il, 1) * qs(il, 1)) |
131 |
else |
132 |
qp(il, i - 1) = qp(il, i) + (cpd * (t(il, i) & |
133 |
- t(il, i - 1)) + gz(il, i) - gz(il, i - 1)) / lv(il, i) |
134 |
qp(il, i - 1) = 0.5 * (qp(il, i - 1) + q(il, i - 1)) |
135 |
qp(il, i - 1) = min(qp(il, i - 1), qs(il, i - 1)) |
136 |
qp(il, i - 1) = max(qp(il, i - 1), 0.) |
137 |
afac1 = p(il, i) * (qs(il, i) - qp(il, i)) & |
138 |
/ (1e4 + 2000. * p(il, i) * qs(il, i)) |
139 |
afac2 = p(il, i - 1) * (qs(il, i - 1) - qp(il, i - 1)) & |
140 |
/ (1e4 + 2000. * p(il, i - 1) * qs(il, i - 1)) |
141 |
afac = 0.5 * (afac1 + afac2) |
142 |
endif |
143 |
|
144 |
if (i == inb(il)) afac = 0. |
145 |
afac = max(afac, 0.) |
146 |
bfac = 1. / (sigd * wt(il, i)) |
147 |
|
148 |
! prise en compte de la variation progressive de sigt dans |
149 |
! les couches icb et icb - 1: |
150 |
! pour plcl < ph(i + 1), pr1 = 0 & pr2 = 1 |
151 |
! pour plcl > ph(i), pr1 = 1 & pr2 = 0 |
152 |
! pour ph(i + 1) < plcl < ph(i), pr1 est la proportion a cheval |
153 |
! sur le nuage, et pr2 est la proportion sous la base du |
154 |
! nuage. |
155 |
pr1 = (plcl(il) - ph(il, i + 1)) / (ph(il, i) - ph(il, i + 1)) |
156 |
pr1 = max(0., min(1., pr1)) |
157 |
pr2 = (ph(il, i) - plcl(il)) / (ph(il, i) - ph(il, i + 1)) |
158 |
pr2 = max(0., min(1., pr2)) |
159 |
sigt = sigp(il, i) * pr1 + pr2 |
160 |
|
161 |
b6 = bfac * 50. * sigd * (ph(il, i) - ph(il, i + 1)) * sigt & |
162 |
* afac |
163 |
c6 = water(il, i + 1) + bfac * wdtrain(il) - 50. * sigd * bfac & |
164 |
* (ph(il, i) - ph(il, i + 1)) * evap(il, i + 1) |
165 |
if (c6 > 0.) then |
166 |
revap = 0.5 * (- b6 + sqrt(b6 * b6 + 4. * c6)) |
167 |
evap(il, i) = sigt * afac * revap |
168 |
water(il, i) = revap * revap |
169 |
else |
170 |
evap(il, i) = - evap(il, i + 1) + 0.02 * (wdtrain(il) & |
171 |
+ sigd * wt(il, i) * water(il, i + 1)) & |
172 |
/ (sigd * (ph(il, i) - ph(il, i + 1))) |
173 |
end if |
174 |
|
175 |
! calculate precipitating downdraft mass flux under |
176 |
! hydrostatic approximation |
177 |
|
178 |
if (i /= 1) then |
179 |
tevap = max(0., evap(il, i)) |
180 |
delth = max(0.001, (th(il, i) - th(il, i - 1))) |
181 |
mp(il, i) = 100. * ginv * lvcp(il, i) * sigd * tevap & |
182 |
* (p(il, i - 1) - p(il, i)) / delth |
183 |
|
184 |
! If hydrostatic assumption fails, solve cubic |
185 |
! difference equation for downdraft theta and mass |
186 |
! flux from two simultaneous differential equations |
187 |
|
188 |
amfac = sigd * sigd * 70. * ph(il, i) & |
189 |
* (p(il, i - 1) - p(il, i)) & |
190 |
* (th(il, i) - th(il, i - 1)) / (tv(il, i) * th(il, i)) |
191 |
amp2 = abs(mp(il, i + 1) * mp(il, i + 1) - mp(il, i) & |
192 |
* mp(il, i)) |
193 |
|
194 |
if (amp2 > 0.1 * amfac) then |
195 |
xf = 100. * sigd * sigd * sigd * (ph(il, i) & |
196 |
- ph(il, i + 1)) |
197 |
tf = b(il, i) - 5. * (th(il, i) - th(il, i - 1)) & |
198 |
* t(il, i) / (lvcp(il, i) * sigd * th(il, i)) |
199 |
af = xf * tf + mp(il, i + 1) * mp(il, i + 1) * tinv |
200 |
bf = 2. * (tinv * mp(il, i + 1))**3 + tinv & |
201 |
* mp(il, i + 1) * xf * tf + 50. * (p(il, i - 1) & |
202 |
- p(il, i)) * xf * tevap |
203 |
fac2 = 1. |
204 |
if (bf < 0.) fac2 = - 1. |
205 |
bf = abs(bf) |
206 |
ur = 0.25 * bf * bf - af * af * af * tinv * tinv * tinv |
207 |
|
208 |
if (ur >= 0.) then |
209 |
sru = sqrt(ur) |
210 |
fac = 1. |
211 |
if ((0.5 * bf - sru) < 0.) fac = - 1. |
212 |
mp(il, i) = mp(il, i + 1) * tinv & |
213 |
+ (0.5 * bf + sru)**tinv & |
214 |
+ fac * (abs(0.5 * bf - sru))**tinv |
215 |
else |
216 |
d = atan(2. * sqrt(- ur) / (bf + 1e-28)) |
217 |
if (fac2 < 0.)d = 3.14159 - d |
218 |
mp(il, i) = mp(il, i + 1) * tinv + 2. & |
219 |
* sqrt(af * tinv) * cos(d * tinv) |
220 |
endif |
221 |
|
222 |
mp(il, i) = max(0., mp(il, i)) |
223 |
|
224 |
! Il y a vraisemblablement une erreur dans la |
225 |
! ligne suivante : il faut diviser par (mp(il, |
226 |
! i) * sigd * grav) et non par (mp(il, i) + sigd |
227 |
! * 0.1). Et il faut bien revoir les facteurs |
228 |
! 100. |
229 |
b(il, i - 1) = b(il, i) + 100. * (p(il, i - 1) & |
230 |
- p(il, i)) * tevap / (mp(il, i) + sigd * 0.1) & |
231 |
- 10. * (th(il, i) - th(il, i - 1)) * t(il, i) & |
232 |
/ (lvcp(il, i) * sigd * th(il, i)) |
233 |
b(il, i - 1) = max(b(il, i - 1), 0.) |
234 |
endif |
235 |
|
236 |
! limit magnitude of mp(i) to meet cfl condition |
237 |
ampmax = 2. * (ph(il, i) - ph(il, i + 1)) * delti |
238 |
amp2 = 2. * (ph(il, i - 1) - ph(il, i)) * delti |
239 |
ampmax = min(ampmax, amp2) |
240 |
mp(il, i) = min(mp(il, i), ampmax) |
241 |
|
242 |
! force mp to decrease linearly to zero |
243 |
! between cloud base and the surface |
244 |
if (p(il, i) > p(il, icb(il))) mp(il, i) = mp(il, icb(il)) & |
245 |
* (p(il, 1) - p(il, i)) / (p(il, 1) - p(il, icb(il))) |
246 |
endif ! i == 1 |
247 |
|
248 |
! find mixing ratio of precipitating downdraft |
249 |
|
250 |
if (i /= inb(il)) then |
251 |
qp(il, i) = q(il, i) |
252 |
|
253 |
if (mp(il, i) > mp(il, i + 1)) then |
254 |
qp(il, i) = qp(il, i + 1) * mp(il, i + 1) + q(il, i) & |
255 |
* (mp(il, i) - mp(il, i + 1)) + 100. * ginv & |
256 |
* 0.5 * sigd * (ph(il, i) - ph(il, i + 1)) & |
257 |
* (evap(il, i + 1) + evap(il, i)) |
258 |
qp(il, i) = qp(il, i) / mp(il, i) |
259 |
up(il, i) = up(il, i + 1) * mp(il, i + 1) + u(il, i) & |
260 |
* (mp(il, i) - mp(il, i + 1)) |
261 |
up(il, i) = up(il, i) / mp(il, i) |
262 |
vp(il, i) = vp(il, i + 1) * mp(il, i + 1) + v(il, i) & |
263 |
* (mp(il, i) - mp(il, i + 1)) |
264 |
vp(il, i) = vp(il, i) / mp(il, i) |
265 |
else |
266 |
if (mp(il, i + 1) > 1e-16) then |
267 |
qp(il, i) = qp(il, i + 1) + 100. * ginv * 0.5 * sigd & |
268 |
* (ph(il, i) - ph(il, i + 1)) & |
269 |
* (evap(il, i + 1) + evap(il, i)) / mp(il, i + 1) |
270 |
up(il, i) = up(il, i + 1) |
271 |
vp(il, i) = vp(il, i + 1) |
272 |
endif |
273 |
endif |
274 |
|
275 |
qp(il, i) = min(qp(il, i), qs(il, i)) |
276 |
qp(il, i) = max(qp(il, i), 0.) |
277 |
endif |
278 |
endif |
279 |
end do |
280 |
end if |
281 |
end DO downdraft_loop |
282 |
|
283 |
end SUBROUTINE cv30_unsat |
284 |
|
285 |
end module cv30_unsat_m |