/[lmdze]/trunk/Sources/phylmd/CV30_routines/cv30_yield.f
ViewVC logotype

Annotation of /trunk/Sources/phylmd/CV30_routines/cv30_yield.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 198 - (hide annotations)
Tue May 31 16:17:35 2016 UTC (7 years, 11 months ago) by guez
File size: 19673 byte(s)
Removed variables nk1 and nk in cv_driver and below. These arrays were
just equal to the constant minorig. (This is also the case in LMDZ.)

In cv_thermo, removed some variables which were copies of variables of
suphec_m. Changed some variables to named constants.

1 guez 185 module cv30_yield_m
2 guez 47
3 guez 97 implicit none
4 guez 47
5 guez 97 contains
6 guez 47
7 guez 188 SUBROUTINE cv30_yield(icb, inb, delt, t, rr, u, v, gz, p, ph, h, hp, lv, &
8     cpn, th, ep, clw, m, tp, mp, rp, up, vp, wt, water, evap, b, ment, &
9     qent, uent, vent, nent, elij, sig, tv, tvp, iflag, precip, VPrecip, &
10 guez 189 ft, fr, fu, fv, upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc)
11 guez 47
12 guez 195 ! Tendencies, precipitation, variables of interface with other
13     ! processes, etc.
14    
15 guez 187 use conema3_m, only: iflag_clw
16 guez 195 use cv30_param_m, only: minorig, nl, sigd
17 guez 198 use cv_thermo_m, only: cl, cpd, cpv, rowl, rrd, rrv
18 guez 188 USE dimphy, ONLY: klev, klon
19 guez 198 use SUPHEC_M, only: rg
20 guez 187
21 guez 97 ! inputs:
22 guez 188 integer, intent(in):: icb(:), inb(:) ! (ncum)
23 guez 97 real, intent(in):: delt
24 guez 188 real t(klon, klev), rr(klon, klev), u(klon, klev), v(klon, klev)
25     real gz(klon, klev)
26     real p(klon, klev)
27     real ph(klon, klev + 1), h(klon, klev), hp(klon, klev)
28     real lv(klon, klev), cpn(klon, klev)
29     real th(klon, klev)
30     real ep(klon, klev), clw(klon, klev)
31     real m(klon, klev)
32     real tp(klon, klev)
33     real mp(klon, klev), rp(klon, klev), up(klon, klev)
34 guez 195 real, intent(in):: vp(:, 2:) ! (ncum, 2:nl)
35     real, intent(in):: wt(:, :) ! (ncum, nl - 1)
36 guez 189 real, intent(in):: water(:, :), evap(:, :) ! (ncum, nl)
37     real, intent(in):: b(:, :) ! (ncum, nl - 1)
38 guez 188 real ment(klon, klev, klev), qent(klon, klev, klev), uent(klon, klev, klev)
39     real vent(klon, klev, klev)
40     integer nent(klon, klev)
41     real elij(klon, klev, klev)
42     real sig(klon, klev)
43     real tv(klon, klev), tvp(klon, klev)
44 guez 47
45 guez 196 integer, intent(out):: iflag(:) ! (ncum)
46 guez 47
47 guez 97 ! outputs:
48 guez 188 real precip(klon)
49     real VPrecip(klon, klev + 1)
50     real ft(klon, klev), fr(klon, klev), fu(klon, klev), fv(klon, klev)
51     real upwd(klon, klev), dnwd(klon, klev)
52     real dnwd0(klon, klev)
53     real ma(klon, klev)
54     real mike(klon, klev)
55     real tls(klon, klev), tps(klon, klev)
56     real qcondc(klon, klev)
57 guez 97
58 guez 188 ! Local:
59 guez 195 real, parameter:: delta = 0.01 ! interface cloud parameterization
60 guez 188 integer ncum
61 guez 187 integer i, k, il, n, j, num1
62 guez 97 real rat, awat, delti
63 guez 105 real ax, bx, cx, dx
64 guez 97 real cpinv, rdcp, dpinv
65 guez 188 real lvcp(klon, klev)
66     real am(klon), work(klon), ad(klon), amp1(klon)
67     real up1(klon, klev, klev), dn1(klon, klev, klev)
68     real asum(klon), bsum(klon), csum(klon), dsum(klon)
69     real qcond(klon, klev), nqcond(klon, klev), wa(klon, klev)
70     real siga(klon, klev), sax(klon, klev), mac(klon, klev)
71 guez 47
72 guez 97 !-------------------------------------------------------------
73 guez 47
74 guez 188 ncum = size(icb)
75 guez 196 iflag = 0
76 guez 188
77 guez 97 ! initialization:
78 guez 47
79 guez 188 delti = 1.0 / delt
80 guez 47
81 guez 188 do il = 1, ncum
82     precip(il) = 0.0
83     VPrecip(il, klev + 1) = 0.
84 guez 97 enddo
85 guez 47
86 guez 188 do i = 1, klev
87     do il = 1, ncum
88     VPrecip(il, i) = 0.0
89     ft(il, i) = 0.0
90     fr(il, i) = 0.0
91     fu(il, i) = 0.0
92     fv(il, i) = 0.0
93     qcondc(il, i) = 0.0
94     qcond(il, i) = 0.0
95     nqcond(il, i) = 0.0
96 guez 47 enddo
97 guez 97 enddo
98 guez 47
99 guez 188 do i = 1, nl
100     do il = 1, ncum
101     lvcp(il, i) = lv(il, i) / cpn(il, i)
102 guez 47 enddo
103 guez 97 enddo
104 guez 47
105 guez 188 ! calculate surface precipitation in mm / day
106 guez 47
107 guez 188 do il = 1, ncum
108     if (ep(il, inb(il)) >= 1e-4) precip(il) = wt(il, 1) * sigd &
109 guez 198 * water(il, 1) * 86400. * 1000. / (rowl * rg)
110 guez 97 enddo
111 guez 47
112 guez 188 ! CALCULATE VERTICAL PROFILE OF PRECIPITATIONs IN kg / m2 / s ===
113 guez 187
114 guez 97 ! MAF rajout pour lessivage
115 guez 188 do k = 1, nl - 1
116     do il = 1, ncum
117     if (k <= inb(il)) VPrecip(il, k) = wt(il, k) * sigd * water(il, k) &
118 guez 198 / rg
119 guez 47 end do
120 guez 97 end do
121 guez 187
122 guez 188 ! calculate tendencies of lowest level potential temperature
123     ! and mixing ratio
124 guez 187
125 guez 188 do il = 1, ncum
126     work(il) = 1.0 / (ph(il, 1) - ph(il, 2))
127     am(il) = 0.0
128 guez 97 enddo
129 guez 47
130 guez 188 do k = 2, nl
131     do il = 1, ncum
132     if (k <= inb(il)) am(il) = am(il) + m(il, k)
133 guez 47 enddo
134 guez 97 enddo
135 guez 47
136 guez 188 do il = 1, ncum
137     ! Consist vect:
138 guez 198 if (0.01 * rg * work(il) * am(il) >= delti) iflag(il) = 1
139 guez 47
140 guez 198 ft(il, 1) = 0.01 * rg * work(il) * am(il) * (t(il, 2) - t(il, 1) &
141 guez 188 + (gz(il, 2) - gz(il, 1)) / cpn(il, 1))
142 guez 47
143 guez 188 ft(il, 1) = ft(il, 1) - 0.5 * lvcp(il, 1) * sigd * (evap(il, 1) &
144     + evap(il, 2))
145 guez 47
146 guez 198 ft(il, 1) = ft(il, 1) - 0.009 * rg * sigd * mp(il, 2) &
147 guez 188 * t(il, 1) * b(il, 1) * work(il)
148 guez 47
149 guez 188 ft(il, 1) = ft(il, 1) + 0.01 * sigd * wt(il, 1) * (cl - cpd) &
150     * water(il, 2) * (t(il, 2) - t(il, 1)) * work(il) / cpn(il, 1)
151 guez 47
152 guez 187 !jyg1 Correction pour mieux conserver l'eau (conformite avec CONVECT4.3)
153 guez 198 ! (sb: pour l'instant, on ne fait que le chgt concernant rg, pas evap)
154     fr(il, 1) = 0.01 * rg * mp(il, 2) * (rp(il, 2) - rr(il, 1)) &
155 guez 188 * work(il) + sigd * 0.5 * (evap(il, 1) + evap(il, 2))
156     ! + tard : + sigd * evap(il, 1)
157 guez 47
158 guez 198 fr(il, 1) = fr(il, 1) + 0.01 * rg * am(il) * (rr(il, 2) - rr(il, 1)) &
159 guez 188 * work(il)
160 guez 47
161 guez 198 fu(il, 1) = fu(il, 1) + 0.01 * rg * work(il) * (mp(il, 2) &
162 guez 188 * (up(il, 2) - u(il, 1)) + am(il) * (u(il, 2) - u(il, 1)))
163 guez 198 fv(il, 1) = fv(il, 1) + 0.01 * rg * work(il) * (mp(il, 2) &
164 guez 188 * (vp(il, 2) - v(il, 1)) + am(il) * (v(il, 2) - v(il, 1)))
165 guez 97 enddo ! il
166 guez 47
167 guez 188 do j = 2, nl
168     do il = 1, ncum
169 guez 187 if (j <= inb(il)) then
170 guez 198 fr(il, 1) = fr(il, 1) + 0.01 * rg * work(il) * ment(il, j, 1) &
171 guez 188 * (qent(il, j, 1) - rr(il, 1))
172 guez 198 fu(il, 1) = fu(il, 1) + 0.01 * rg * work(il) * ment(il, j, 1) &
173 guez 188 * (uent(il, j, 1) - u(il, 1))
174 guez 198 fv(il, 1) = fv(il, 1) + 0.01 * rg * work(il) * ment(il, j, 1) &
175 guez 188 * (vent(il, j, 1) - v(il, 1))
176     endif
177 guez 47 enddo
178 guez 97 enddo
179 guez 47
180 guez 188 ! calculate tendencies of potential temperature and mixing ratio
181     ! at levels above the lowest level
182 guez 47
183 guez 188 ! first find the net saturated updraft and downdraft mass fluxes
184     ! through each level
185 guez 47
186 guez 188 loop_i: do i = 2, nl - 1
187     num1 = 0
188 guez 187
189 guez 188 do il = 1, ncum
190     if (i <= inb(il)) num1 = num1 + 1
191 guez 47 enddo
192    
193 guez 188 if (num1 > 0) then
194     amp1(:ncum) = 0.
195     ad(:ncum) = 0.
196 guez 47
197 guez 188 do k = i + 1, nl + 1
198     do il = 1, ncum
199     if (i <= inb(il) .and. k <= (inb(il) + 1)) then
200     amp1(il) = amp1(il) + m(il, k)
201     endif
202     end do
203 guez 97 end do
204 guez 47
205 guez 188 do k = 1, i
206     do j = i + 1, nl + 1
207     do il = 1, ncum
208     if (i <= inb(il) .and. j <= (inb(il) + 1)) then
209     amp1(il) = amp1(il) + ment(il, k, j)
210     endif
211     end do
212 guez 97 end do
213     end do
214 guez 47
215 guez 188 do k = 1, i - 1
216     do j = i, nl + 1 ! newvecto: nl au lieu nl + 1?
217     do il = 1, ncum
218     if (i <= inb(il) .and. j <= inb(il)) then
219     ad(il) = ad(il) + ment(il, j, k)
220     endif
221     end do
222 guez 97 end do
223     end do
224 guez 47
225 guez 188 do il = 1, ncum
226     if (i <= inb(il)) then
227     dpinv = 1.0 / (ph(il, i) - ph(il, i + 1))
228     cpinv = 1.0 / cpn(il, i)
229 guez 47
230 guez 188 ! Vecto:
231 guez 198 if (0.01 * rg * dpinv * amp1(il) >= delti) iflag(il) = 1
232 guez 47
233 guez 198 ft(il, i) = 0.01 * rg * dpinv * (amp1(il) * (t(il, i + 1) &
234 guez 188 - t(il, i) + (gz(il, i + 1) - gz(il, i)) * cpinv) &
235     - ad(il) * (t(il, i) - t(il, i - 1) + (gz(il, i) &
236     - gz(il, i - 1)) * cpinv)) - 0.5 * sigd * lvcp(il, i) &
237     * (evap(il, i) + evap(il, i + 1))
238     rat = cpn(il, i - 1) * cpinv
239 guez 198 ft(il, i) = ft(il, i) - 0.009 * rg * sigd * (mp(il, i + 1) &
240 guez 188 * t(il, i) * b(il, i) - mp(il, i) * t(il, i - 1) * rat &
241     * b(il, i - 1)) * dpinv
242 guez 198 ft(il, i) = ft(il, i) + 0.01 * rg * dpinv * ment(il, i, i) &
243 guez 188 * (hp(il, i) - h(il, i) + t(il, i) * (cpv - cpd) &
244     * (rr(il, i) - qent(il, i, i))) * cpinv
245 guez 47
246 guez 188 ft(il, i) = ft(il, i) + 0.01 * sigd * wt(il, i) * (cl - cpd) &
247     * water(il, i + 1) * (t(il, i + 1) - t(il, i)) * dpinv &
248     * cpinv
249 guez 47
250 guez 198 fr(il, i) = 0.01 * rg * dpinv * (amp1(il) * (rr(il, i + 1) &
251 guez 188 - rr(il, i)) - ad(il) * (rr(il, i) - rr(il, i - 1)))
252 guez 198 fu(il, i) = fu(il, i) + 0.01 * rg * dpinv * (amp1(il) &
253 guez 188 * (u(il, i + 1) - u(il, i)) - ad(il) * (u(il, i) &
254     - u(il, i - 1)))
255 guez 198 fv(il, i) = fv(il, i) + 0.01 * rg * dpinv * (amp1(il) &
256 guez 188 * (v(il, i + 1) - v(il, i)) - ad(il) * (v(il, i) &
257     - v(il, i - 1)))
258     endif
259     end do
260 guez 47
261 guez 188 do k = 1, i - 1
262     do il = 1, ncum
263     if (i <= inb(il)) then
264     dpinv = 1.0 / (ph(il, i) - ph(il, i + 1))
265     cpinv = 1.0 / cpn(il, i)
266 guez 47
267 guez 188 awat = elij(il, k, i) - (1. - ep(il, i)) * clw(il, i)
268     awat = amax1(awat, 0.0)
269 guez 47
270 guez 198 fr(il, i) = fr(il, i) + 0.01 * rg * dpinv &
271 guez 188 * ment(il, k, i) * (qent(il, k, i) - awat - rr(il, i))
272 guez 198 fu(il, i) = fu(il, i) + 0.01 * rg * dpinv &
273 guez 188 * ment(il, k, i) * (uent(il, k, i) - u(il, i))
274 guez 198 fv(il, i) = fv(il, i) + 0.01 * rg * dpinv &
275 guez 188 * ment(il, k, i) * (vent(il, k, i) - v(il, i))
276 guez 47
277 guez 188 ! (saturated updrafts resulting from mixing)
278     qcond(il, i) = qcond(il, i) + (elij(il, k, i) - awat)
279     nqcond(il, i) = nqcond(il, i) + 1.
280     endif ! i
281     end do
282 guez 97 end do
283 guez 47
284 guez 188 do k = i, nl + 1
285     do il = 1, ncum
286     if (i <= inb(il) .and. k <= inb(il)) then
287     dpinv = 1.0 / (ph(il, i) - ph(il, i + 1))
288     cpinv = 1.0 / cpn(il, i)
289 guez 47
290 guez 198 fr(il, i) = fr(il, i) + 0.01 * rg * dpinv &
291 guez 188 * ment(il, k, i) * (qent(il, k, i) - rr(il, i))
292 guez 198 fu(il, i) = fu(il, i) + 0.01 * rg * dpinv &
293 guez 188 * ment(il, k, i) * (uent(il, k, i) - u(il, i))
294 guez 198 fv(il, i) = fv(il, i) + 0.01 * rg * dpinv &
295 guez 188 * ment(il, k, i) * (vent(il, k, i) - v(il, i))
296     endif
297     end do
298 guez 97 end do
299 guez 47
300 guez 188 do il = 1, ncum
301     if (i <= inb(il)) then
302     dpinv = 1.0 / (ph(il, i) - ph(il, i + 1))
303     cpinv = 1.0 / cpn(il, i)
304 guez 47
305 guez 188 ! sb: on ne fait pas encore la correction permettant de mieux
306     ! conserver l'eau:
307     fr(il, i) = fr(il, i) + 0.5 * sigd * (evap(il, i) &
308 guez 198 + evap(il, i + 1)) + 0.01 * rg * (mp(il, i + 1) &
309 guez 188 * (rp(il, i + 1) - rr(il, i)) - mp(il, i) * (rp(il, i) &
310     - rr(il, i - 1))) * dpinv
311 guez 47
312 guez 198 fu(il, i) = fu(il, i) + 0.01 * rg * (mp(il, i + 1) &
313 guez 188 * (up(il, i + 1) - u(il, i)) - mp(il, i) * (up(il, i) &
314     - u(il, i - 1))) * dpinv
315 guez 198 fv(il, i) = fv(il, i) + 0.01 * rg * (mp(il, i + 1) &
316 guez 188 * (vp(il, i + 1) - v(il, i)) - mp(il, i) * (vp(il, i) &
317     - v(il, i - 1))) * dpinv
318     endif
319     end do
320 guez 47
321 guez 188 ! sb: interface with the cloud parameterization:
322 guez 47
323 guez 188 do k = i + 1, nl
324     do il = 1, ncum
325     if (k <= inb(il) .and. i <= inb(il)) then
326     ! (saturated downdrafts resulting from mixing)
327     qcond(il, i) = qcond(il, i) + elij(il, k, i)
328     nqcond(il, i) = nqcond(il, i) + 1.
329     endif
330     enddo
331     enddo
332 guez 47
333 guez 188 ! (particular case: no detraining level is found)
334     do il = 1, ncum
335     if (i <= inb(il) .and. nent(il, i) == 0) then
336     qcond(il, i) = qcond(il, i) + (1. - ep(il, i)) * clw(il, i)
337     nqcond(il, i) = nqcond(il, i) + 1.
338     endif
339     enddo
340 guez 47
341 guez 188 do il = 1, ncum
342     if (i <= inb(il) .and. nqcond(il, i) /= 0.) then
343     qcond(il, i) = qcond(il, i) / nqcond(il, i)
344     endif
345     enddo
346     end if
347     end do loop_i
348 guez 47
349 guez 188 ! move the detrainment at level inb down to level inb - 1
350     ! in such a way as to preserve the vertically
351     ! integrated enthalpy and water tendencies
352 guez 47
353 guez 188 do il = 1, ncum
354     ax = 0.1 * ment(il, inb(il), inb(il)) * (hp(il, inb(il)) &
355     - h(il, inb(il)) + t(il, inb(il)) * (cpv - cpd) &
356     * (rr(il, inb(il)) - qent(il, inb(il), inb(il)))) &
357     / (cpn(il, inb(il)) * (ph(il, inb(il)) - ph(il, inb(il) + 1)))
358     ft(il, inb(il)) = ft(il, inb(il)) - ax
359     ft(il, inb(il) - 1) = ft(il, inb(il) - 1) + ax * cpn(il, inb(il)) &
360     * (ph(il, inb(il)) - ph(il, inb(il) + 1)) / (cpn(il, inb(il) - 1) &
361     * (ph(il, inb(il) - 1) - ph(il, inb(il))))
362 guez 47
363 guez 188 bx = 0.1 * ment(il, inb(il), inb(il)) * (qent(il, inb(il), inb(il)) &
364     - rr(il, inb(il))) / (ph(il, inb(il)) - ph(il, inb(il) + 1))
365     fr(il, inb(il)) = fr(il, inb(il)) - bx
366     fr(il, inb(il) - 1) = fr(il, inb(il) - 1) &
367     + bx * (ph(il, inb(il)) - ph(il, inb(il) + 1)) &
368     / (ph(il, inb(il) - 1) - ph(il, inb(il)))
369 guez 47
370 guez 188 cx = 0.1 * ment(il, inb(il), inb(il)) * (uent(il, inb(il), inb(il)) &
371     - u(il, inb(il))) / (ph(il, inb(il)) - ph(il, inb(il) + 1))
372     fu(il, inb(il)) = fu(il, inb(il)) - cx
373     fu(il, inb(il) - 1) = fu(il, inb(il) - 1) &
374     + cx * (ph(il, inb(il)) - ph(il, inb(il) + 1)) &
375     / (ph(il, inb(il) - 1) - ph(il, inb(il)))
376 guez 47
377 guez 188 dx = 0.1 * ment(il, inb(il), inb(il)) * (vent(il, inb(il), inb(il)) &
378     - v(il, inb(il))) / (ph(il, inb(il)) - ph(il, inb(il) + 1))
379     fv(il, inb(il)) = fv(il, inb(il)) - dx
380     fv(il, inb(il) - 1) = fv(il, inb(il) - 1) &
381     + dx * (ph(il, inb(il)) - ph(il, inb(il) + 1)) &
382     / (ph(il, inb(il) - 1) - ph(il, inb(il)))
383 guez 47
384 guez 97 end do
385 guez 47
386 guez 188 ! homoginize tendencies below cloud base
387 guez 187
388 guez 188 do il = 1, ncum
389     asum(il) = 0.0
390     bsum(il) = 0.0
391     csum(il) = 0.0
392     dsum(il) = 0.0
393 guez 97 enddo
394 guez 47
395 guez 188 do i = 1, nl
396     do il = 1, ncum
397     if (i <= (icb(il) - 1)) then
398     asum(il) = asum(il) + ft(il, i) * (ph(il, i) - ph(il, i + 1))
399     bsum(il) = bsum(il) + fr(il, i) * (lv(il, i) + (cl - cpd) &
400     * (t(il, i) - t(il, 1))) * (ph(il, i) - ph(il, i + 1))
401     csum(il) = csum(il) + (lv(il, i) + (cl - cpd) * (t(il, i) &
402     - t(il, 1))) * (ph(il, i) - ph(il, i + 1))
403     dsum(il) = dsum(il) + t(il, i) * (ph(il, i) - ph(il, i + 1)) &
404     / th(il, i)
405 guez 97 endif
406 guez 47 enddo
407 guez 97 enddo
408 guez 47
409 guez 188 do i = 1, nl
410     do il = 1, ncum
411     if (i <= (icb(il) - 1)) then
412     ft(il, i) = asum(il) * t(il, i) / (th(il, i) * dsum(il))
413     fr(il, i) = bsum(il) / csum(il)
414 guez 97 endif
415 guez 47 enddo
416 guez 97 enddo
417 guez 47
418 guez 188 ! reset counter and return
419 guez 187
420 guez 188 do il = 1, ncum
421     sig(il, klev) = 2.0
422 guez 97 enddo
423 guez 47
424 guez 188 do i = 1, klev
425     do il = 1, ncum
426     upwd(il, i) = 0.0
427     dnwd(il, i) = 0.0
428 guez 47 enddo
429 guez 97 enddo
430 guez 47
431 guez 188 do i = 1, nl
432     do il = 1, ncum
433     dnwd0(il, i) = - mp(il, i)
434 guez 47 enddo
435 guez 97 enddo
436 guez 188 do i = nl + 1, klev
437     do il = 1, ncum
438     dnwd0(il, i) = 0.
439 guez 47 enddo
440 guez 97 enddo
441 guez 47
442 guez 188 do i = 1, nl
443     do il = 1, ncum
444 guez 187 if (i >= icb(il) .and. i <= inb(il)) then
445 guez 188 upwd(il, i) = 0.0
446     dnwd(il, i) = 0.0
447 guez 97 endif
448 guez 47 enddo
449 guez 97 enddo
450 guez 47
451 guez 188 do i = 1, nl
452     do k = 1, nl
453     do il = 1, ncum
454     up1(il, k, i) = 0.0
455     dn1(il, k, i) = 0.0
456 guez 97 enddo
457 guez 47 enddo
458 guez 97 enddo
459 guez 47
460 guez 188 do i = 1, nl
461     do k = i, nl
462     do n = 1, i - 1
463     do il = 1, ncum
464 guez 187 if (i >= icb(il).and.i <= inb(il).and.k <= inb(il)) then
465 guez 188 up1(il, k, i) = up1(il, k, i) + ment(il, n, k)
466     dn1(il, k, i) = dn1(il, k, i) - ment(il, k, n)
467 guez 97 endif
468     enddo
469     enddo
470 guez 47 enddo
471 guez 97 enddo
472 guez 47
473 guez 188 do i = 2, nl
474     do k = i, nl
475     do il = 1, ncum
476 guez 187 if (i <= inb(il).and.k <= inb(il)) then
477 guez 188 upwd(il, i) = upwd(il, i) + m(il, k) + up1(il, k, i)
478     dnwd(il, i) = dnwd(il, i) + dn1(il, k, i)
479 guez 97 endif
480     enddo
481 guez 47 enddo
482 guez 97 enddo
483 guez 47
484 guez 97 !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
485 guez 187 ! determination de la variation de flux ascendant entre
486     ! deux niveau non dilue mike
487 guez 97 !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
488 guez 47
489 guez 188 do i = 1, nl
490     do il = 1, ncum
491     mike(il, i) = m(il, i)
492 guez 47 enddo
493 guez 97 enddo
494 guez 47
495 guez 188 do i = nl + 1, klev
496     do il = 1, ncum
497     mike(il, i) = 0.
498 guez 47 enddo
499 guez 97 enddo
500 guez 47
501 guez 188 do i = 1, klev
502     do il = 1, ncum
503     ma(il, i) = 0
504 guez 47 enddo
505 guez 97 enddo
506 guez 47
507 guez 188 do i = 1, nl
508     do j = i, nl
509     do il = 1, ncum
510     ma(il, i) = ma(il, i) + m(il, j)
511 guez 97 enddo
512 guez 47 enddo
513 guez 97 enddo
514 guez 47
515 guez 188 do i = nl + 1, klev
516     do il = 1, ncum
517     ma(il, i) = 0.
518 guez 47 enddo
519 guez 97 enddo
520 guez 47
521 guez 188 do i = 1, nl
522     do il = 1, ncum
523     if (i <= (icb(il) - 1)) then
524     ma(il, i) = 0
525 guez 97 endif
526 guez 47 enddo
527 guez 97 enddo
528 guez 47
529 guez 97 !cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
530 guez 187 ! icb represente de niveau ou se trouve la
531     ! base du nuage, et inb le top du nuage
532 guez 97 !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
533 guez 47
534 guez 188 do i = 1, klev
535     DO il = 1, ncum
536     rdcp = (rrd * (1. - rr(il, i)) - rr(il, i) * rrv) &
537     / (cpd * (1. - rr(il, i)) + rr(il, i) * cpv)
538     tls(il, i) = t(il, i) * (1000.0 / p(il, i))**rdcp
539     tps(il, i) = tp(il, i)
540 guez 97 end DO
541     enddo
542 guez 47
543 guez 188 ! diagnose the in-cloud mixing ratio
544     ! of condensed water
545     !
546 guez 47
547 guez 188 do i = 1, klev
548     do il = 1, ncum
549     mac(il, i) = 0.0
550     wa(il, i) = 0.0
551     siga(il, i) = 0.0
552     sax(il, i) = 0.0
553     enddo
554     enddo
555 guez 47
556 guez 188 do i = minorig, nl
557     do k = i + 1, nl + 1
558     do il = 1, ncum
559     if (i <= inb(il) .and. k <= (inb(il) + 1)) then
560     mac(il, i) = mac(il, i) + m(il, k)
561     endif
562     enddo
563     enddo
564     enddo
565 guez 47
566 guez 188 do i = 1, nl
567     do j = 1, i
568     do il = 1, ncum
569     if (i >= icb(il) .and. i <= (inb(il) - 1) &
570     .and. j >= icb(il)) then
571     sax(il, i) = sax(il, i) + rrd * (tvp(il, j) - tv(il, j)) &
572     * (ph(il, j) - ph(il, j + 1)) / p(il, j)
573     endif
574     enddo
575     enddo
576     enddo
577 guez 97
578 guez 188 do i = 1, nl
579     do il = 1, ncum
580     if (i >= icb(il) .and. i <= (inb(il) - 1) &
581     .and. sax(il, i) > 0.0) then
582     wa(il, i) = sqrt(2. * sax(il, i))
583     endif
584     enddo
585     enddo
586 guez 47
587 guez 188 do i = 1, nl
588     do il = 1, ncum
589     if (wa(il, i) > 0.0) siga(il, i) = mac(il, i) / wa(il, i) * rrd &
590     * tvp(il, i) / p(il, i) / 100. / delta
591     siga(il, i) = min(siga(il, i), 1.0)
592    
593 guez 187 if (iflag_clw == 0) then
594 guez 188 qcondc(il, i) = siga(il, i) * clw(il, i) * (1. - ep(il, i)) &
595     + (1. - siga(il, i)) * qcond(il, i)
596 guez 187 else if (iflag_clw == 1) then
597 guez 188 qcondc(il, i) = qcond(il, i)
598 guez 97 endif
599 guez 188 enddo
600     enddo
601 guez 47
602 guez 185 end SUBROUTINE cv30_yield
603 guez 97
604 guez 185 end module cv30_yield_m

  ViewVC Help
Powered by ViewVC 1.1.21