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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 194 - (hide annotations)
Thu May 12 14:35:35 2016 UTC (8 years ago) by guez
File size: 10777 byte(s)
Clarified computation of wdtrain in procedure cv30_unsat: sum
intrinsic instead of loop. Changes results.

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

  ViewVC Help
Powered by ViewVC 1.1.21