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

Contents of /trunk/Sources/phylmd/CV30_routines/cv30_mixing.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 195 - (show annotations)
Wed May 18 17:56:44 2016 UTC (8 years ago) by guez
File size: 10903 byte(s)
In cv30_feed, iflag1 is 0 on entry so we can simplify the test for
iflag1 = 7.

In cv30_feed, for the computation of icb, replaced sequential search
(with a useless end of loop on k) by a call to locate.

In CV30 routines, replaced len, nloc, nd, na by klon or
klev. Philosophy: no more generality than actually necessary.

Converted as many variables as possible to named constants in
cv30_param_m and downgraded pbcrit, ptcrit, dtovsh, dpbase, dttrig,
tau, delta to local objects in procedures. spfac, betad and omtrain
are useless and removed.

Instead of filling the array sigp with the constant spfac in
cv30_undilute2, just made sigp a constant in cv30_unsat.

In cv_driver, define as allocatable variables that are only
used on the range (ncum, nl).

1 module cv30_mixing_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE cv30_mixing(icb, nk, inb, t, rr, rs, u, v, h, lv, hp, ep, clw, &
8 m, sig, ment, qent, uent, vent, nent, sij, elij, ments, qents)
9
10 ! MIXING
11
12 ! a faire:
13 ! - changer rr(il, 1) -> qnk(il)
14 ! - vectorisation de la partie normalisation des flux (do 789)
15
16 use cv30_param_m, only: minorig, nl
17 use cv_thermo_m, only: cpd, cpv, rrv
18 USE dimphy, ONLY: klev, klon
19
20 ! inputs:
21 integer, intent(in):: icb(:), nk(:), inb(:) ! (ncum)
22 real t(klon, klev), rr(klon, klev), rs(klon, klev)
23 real u(klon, klev), v(klon, klev)
24 real h(klon, klev), lv(klon, klev), hp(klon, klev)
25 real ep(klon, klev), clw(klon, klev)
26 real m(klon, klev) ! input of convect3
27 real sig(klon, klev)
28
29 ! outputs:
30 real ment(klon, klev, klev), qent(klon, klev, klev)
31 real uent(klon, klev, klev), vent(klon, klev, klev)
32 integer nent(klon, klev)
33 real sij(klon, klev, klev), elij(klon, klev, klev)
34 real ments(klon, klev, klev), qents(klon, klev, klev)
35
36 ! Local:
37 integer ncum, i, j, k, il, im, jm
38 integer num1, num2
39 real rti, bf2, anum, denom, dei, altem, cwat, stemp, qp
40 real alt, smid, sjmin, sjmax, delp, delm
41 real asij(klon), smax(klon), scrit(klon)
42 real asum(klon, klev), bsum(klon, klev), csum(klon, klev)
43 real wgh
44 real zm(klon, klev)
45 logical lwork(klon)
46
47 !-------------------------------------------------------------------------
48
49 ncum = size(icb)
50
51 ! INITIALIZE VARIOUS ARRAYS USED IN THE COMPUTATIONS
52
53 do j = 1, nl
54 do i = 1, ncum
55 nent(i, j) = 0
56 end do
57 end do
58
59 do j = 1, nl
60 do k = 1, nl
61 do i = 1, ncum
62 qent(i, k, j) = rr(i, j)
63 uent(i, k, j) = u(i, j)
64 vent(i, k, j) = v(i, j)
65 elij(i, k, j) = 0.0
66 end do
67 end do
68 end do
69
70 ment(1:ncum, 1:klev, 1:klev) = 0.0
71 sij(1:ncum, 1:klev, 1:klev) = 0.0
72
73 zm(:, :) = 0.
74
75 ! CALCULATE ENTRAINED AIR MASS FLUX (ment), TOTAL WATER MIXING
76 ! RATIO (QENT), TOTAL CONDENSED WATER (elij), AND MIXING
77 ! FRACTION (sij)
78
79 do i = minorig + 1, nl
80
81 do j = minorig, nl
82 do il = 1, ncum
83 if((i >= icb(il)).and.(i <= inb(il)).and. &
84 (j >= (icb(il) - 1)).and.(j <= inb(il)))then
85
86 rti = rr(il, 1) - ep(il, i) * clw(il, i)
87 bf2 = 1. + lv(il, j) * lv(il, j) * rs(il, j) / (rrv &
88 * t(il, j) * t(il, j) * cpd)
89 anum = h(il, j) - hp(il, i) + (cpv - cpd) * t(il, j) * (rti &
90 - rr(il, j))
91 denom = h(il, i) - hp(il, i) + (cpd - cpv) * (rr(il, i) &
92 - rti) * t(il, j)
93 dei = denom
94 if(abs(dei) < 0.01)dei = 0.01
95 sij(il, i, j) = anum / dei
96 sij(il, i, i) = 1.0
97 altem = sij(il, i, j) * rr(il, i) + (1. - sij(il, i, j)) &
98 * rti - rs(il, j)
99 altem = altem / bf2
100 cwat = clw(il, j) * (1. - ep(il, j))
101 stemp = sij(il, i, j)
102 if((stemp < 0.0.or.stemp > 1.0.or.altem > cwat) &
103 .and.j > i)then
104 anum = anum - lv(il, j) * (rti - rs(il, j) - cwat * bf2)
105 denom = denom + lv(il, j) * (rr(il, i) - rti)
106 if(abs(denom) < 0.01)denom = 0.01
107 sij(il, i, j) = anum / denom
108 altem = sij(il, i, j) * rr(il, i) + (1. - sij(il, i, j)) &
109 * rti - rs(il, j)
110 altem = altem - (bf2 - 1.) * cwat
111 end if
112 if(sij(il, i, j) > 0.0.and.sij(il, i, j) < 0.95)then
113 qent(il, i, j) = sij(il, i, j) * rr(il, i) + (1. &
114 - sij(il, i, j)) * rti
115 uent(il, i, j) = sij(il, i, j) * u(il, i) + (1. &
116 - sij(il, i, j)) * u(il, nk(il))
117 vent(il, i, j) = sij(il, i, j) * v(il, i) + (1. &
118 - sij(il, i, j)) * v(il, nk(il))
119 elij(il, i, j) = altem
120 elij(il, i, j) = amax1(0.0, elij(il, i, j))
121 ment(il, i, j) = m(il, i) / (1. - sij(il, i, j))
122 nent(il, i) = nent(il, i) + 1
123 end if
124 sij(il, i, j) = amax1(0.0, sij(il, i, j))
125 sij(il, i, j) = amin1(1.0, sij(il, i, j))
126 endif ! new
127 end do
128 end do
129
130 ! if no air can entrain at level i assume that updraft detrains
131 ! at that level and calculate detrained air flux and properties
132
133 do il = 1, ncum
134 if ((i >= icb(il)).and.(i <= inb(il)).and.(nent(il, i) == 0)) then
135 ment(il, i, i) = m(il, i)
136 qent(il, i, i) = rr(il, nk(il)) - ep(il, i) * clw(il, i)
137 uent(il, i, i) = u(il, nk(il))
138 vent(il, i, i) = v(il, nk(il))
139 elij(il, i, i) = clw(il, i)
140 sij(il, i, i) = 0.0
141 end if
142 end do
143 end do
144
145 ! NORMALIZE ENTRAINED AIR MASS FLUXES
146 ! TO REPRESENT EQUAL PROBABILITIES OF MIXING
147
148 asum = 0.
149 csum = 0.
150
151 do il = 1, ncum
152 lwork(il) = .FALSE.
153 enddo
154
155 DO i = minorig + 1, nl
156
157 num1 = 0
158 do il = 1, ncum
159 if (i >= icb(il) .and. i <= inb(il)) num1 = num1 + 1
160 enddo
161 if (num1 <= 0) cycle
162
163 do il = 1, ncum
164 if (i >= icb(il) .and. i <= inb(il)) then
165 lwork(il) = (nent(il, i) /= 0)
166 qp = rr(il, 1) - ep(il, i) * clw(il, i)
167 anum = h(il, i) - hp(il, i) - lv(il, i) * (qp - rs(il, i)) &
168 + (cpv - cpd) * t(il, i) * (qp - rr(il, i))
169 denom = h(il, i) - hp(il, i) + lv(il, i) * (rr(il, i) - qp) &
170 + (cpd - cpv) * t(il, i) * (rr(il, i) - qp)
171 if(abs(denom) < 0.01)denom = 0.01
172 scrit(il) = anum / denom
173 alt = qp - rs(il, i) + scrit(il) * (rr(il, i) - qp)
174 if(scrit(il) <= 0.0.or.alt <= 0.0)scrit(il) = 1.0
175 smax(il) = 0.0
176 asij(il) = 0.0
177 endif
178 end do
179
180 do j = nl, minorig, - 1
181
182 num2 = 0
183 do il = 1, ncum
184 if (i >= icb(il) .and. i <= inb(il) .and. &
185 j >= (icb(il) - 1) .and. j <= inb(il) &
186 .and. lwork(il)) num2 = num2 + 1
187 enddo
188 if (num2 <= 0) cycle
189
190 do il = 1, ncum
191 if (i >= icb(il) .and. i <= inb(il) .and. &
192 j >= (icb(il) - 1) .and. j <= inb(il) &
193 .and. lwork(il)) then
194
195 if(sij(il, i, j) > 1.0e-16.and.sij(il, i, j) < 0.95)then
196 wgh = 1.0
197 if(j > i)then
198 sjmax = amax1(sij(il, i, j + 1), smax(il))
199 sjmax = amin1(sjmax, scrit(il))
200 smax(il) = amax1(sij(il, i, j), smax(il))
201 sjmin = amax1(sij(il, i, j - 1), smax(il))
202 sjmin = amin1(sjmin, scrit(il))
203 if(sij(il, i, j) < (smax(il) - 1.0e-16))wgh = 0.0
204 smid = amin1(sij(il, i, j), scrit(il))
205 else
206 sjmax = amax1(sij(il, i, j + 1), scrit(il))
207 smid = amax1(sij(il, i, j), scrit(il))
208 sjmin = 0.0
209 if(j > 1)sjmin = sij(il, i, j - 1)
210 sjmin = amax1(sjmin, scrit(il))
211 endif
212 delp = abs(sjmax - smid)
213 delm = abs(sjmin - smid)
214 asij(il) = asij(il) + wgh * (delp + delm)
215 ment(il, i, j) = ment(il, i, j) * (delp + delm) * wgh
216 endif
217 endif
218 end do
219
220 end do
221
222 do il = 1, ncum
223 if (i >= icb(il).and.i <= inb(il).and.lwork(il)) then
224 asij(il) = amax1(1.0e-16, asij(il))
225 asij(il) = 1.0 / asij(il)
226 asum(il, i) = 0.0
227 bsum(il, i) = 0.0
228 csum(il, i) = 0.0
229 endif
230 enddo
231
232 do j = minorig, nl
233 do il = 1, ncum
234 if (i >= icb(il) .and. i <= inb(il) .and. lwork(il) &
235 .and. j >= (icb(il) - 1) .and. j <= inb(il)) then
236 ment(il, i, j) = ment(il, i, j) * asij(il)
237 endif
238 enddo
239 end do
240
241 do j = minorig, nl
242 do il = 1, ncum
243 if (i >= icb(il) .and. i <= inb(il) .and. lwork(il) &
244 .and. j >= (icb(il) - 1) .and. j <= inb(il)) then
245 asum(il, i) = asum(il, i) + ment(il, i, j)
246 ment(il, i, j) = ment(il, i, j) * sig(il, j)
247 bsum(il, i) = bsum(il, i) + ment(il, i, j)
248 endif
249 enddo
250 end do
251
252 do il = 1, ncum
253 if (i >= icb(il).and.i <= inb(il).and.lwork(il)) then
254 bsum(il, i) = amax1(bsum(il, i), 1.0e-16)
255 bsum(il, i) = 1.0 / bsum(il, i)
256 endif
257 enddo
258
259 do j = minorig, nl
260 do il = 1, ncum
261 if (i >= icb(il) .and. i <= inb(il) .and. lwork(il) &
262 .and. j >= (icb(il) - 1) .and. j <= inb(il)) then
263 ment(il, i, j) = ment(il, i, j) * asum(il, i) * bsum(il, i)
264 endif
265 enddo
266 end do
267
268 do j = minorig, nl
269 do il = 1, ncum
270 if (i >= icb(il) .and. i <= inb(il) .and. lwork(il) &
271 .and. j >= (icb(il) - 1) .and. j <= inb(il)) then
272 csum(il, i) = csum(il, i) + ment(il, i, j)
273 endif
274 enddo
275 end do
276
277 do il = 1, ncum
278 if (i >= icb(il) .and. i <= inb(il) .and. lwork(il) &
279 .and. csum(il, i) < m(il, i)) then
280 nent(il, i) = 0
281 ment(il, i, i) = m(il, i)
282 qent(il, i, i) = rr(il, 1) - ep(il, i) * clw(il, i)
283 uent(il, i, i) = u(il, nk(il))
284 vent(il, i, i) = v(il, nk(il))
285 elij(il, i, i) = clw(il, i)
286 sij(il, i, i) = 0.0
287 endif
288 enddo ! il
289
290 end DO
291
292 ! MAF: renormalisation de MENT
293 do jm = 1, klev
294 do im = 1, klev
295 do il = 1, ncum
296 zm(il, im) = zm(il, im) + (1. - sij(il, im, jm)) * ment(il, im, jm)
297 end do
298 end do
299 end do
300
301 do jm = 1, klev
302 do im = 1, klev
303 do il = 1, ncum
304 if(zm(il, im) /= 0.) then
305 ment(il, im, jm) = ment(il, im, jm) * m(il, im) / zm(il, im)
306 endif
307 end do
308 end do
309 end do
310
311 do jm = 1, klev
312 do im = 1, klev
313 do il = 1, ncum
314 qents(il, im, jm) = qent(il, im, jm)
315 ments(il, im, jm) = ment(il, im, jm)
316 end do
317 enddo
318 enddo
319
320 end SUBROUTINE cv30_mixing
321
322 end module cv30_mixing_m

  ViewVC Help
Powered by ViewVC 1.1.21