/[lmdze]/trunk/phylmd/cv_driver.f
ViewVC logotype

Contents of /trunk/phylmd/cv_driver.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 196 - (show annotations)
Mon May 23 13:50:39 2016 UTC (8 years ago) by guez
Original Path: trunk/Sources/phylmd/cv_driver.f
File size: 9520 byte(s)
Removed argument icbmax of cv30_feed, not used in cv_driver (not used
in LMDZ either).

Clearer to have iflag1 = 0 in cv30_feed than in cv_driver. Clearer to
have iflag1 = 42 in cv30_uncompress than in cv_driver.

Removed argument iflag of cv30_compress. Since there is iflag1 = 42 in
cv30_uncompress, iflag1 that comes out of cv_driver is disconnected
from iflag1 computed by cv30_feed and cv30_trigger.

Removed some references to convect3 and convect4 in comments. This
program derives from the convect3 version, we do not need to know
about other versions.

Bug fix in cv30_undilute1: icbs1 was not made >= 2.

1 module cv_driver_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE cv_driver(t1, q1, qs1, u1, v1, p1, ph1, iflag1, ft1, fq1, fu1, &
8 fv1, precip1, VPrecip1, sig1, w01, icb1, inb1, Ma1, upwd1, dnwd1, &
9 dnwd01, qcondc1, cape1, da1, phi1, mp1)
10
11 ! From LMDZ4/libf/phylmd/cv_driver.F, version 1.3, 2005/04/15 12:36:17
12 ! Main driver for convection
13 ! Author: S. Bony, March 2002
14
15 ! Several modules corresponding to different physical processes
16
17 use comconst, only: dtphys
18 use cv30_closure_m, only: cv30_closure
19 use cv30_compress_m, only: cv30_compress
20 use cv30_feed_m, only: cv30_feed
21 use cv30_mixing_m, only: cv30_mixing
22 use cv30_param_m, only: cv30_param, nl
23 use cv30_prelim_m, only: cv30_prelim
24 use cv30_tracer_m, only: cv30_tracer
25 use cv30_trigger_m, only: cv30_trigger
26 use cv30_uncompress_m, only: cv30_uncompress
27 use cv30_undilute1_m, only: cv30_undilute1
28 use cv30_undilute2_m, only: cv30_undilute2
29 use cv30_unsat_m, only: cv30_unsat
30 use cv30_yield_m, only: cv30_yield
31 use cv_thermo_m, only: cv_thermo
32 USE dimphy, ONLY: klev, klon
33
34 real, intent(in):: t1(klon, klev) ! temperature (K)
35 real, intent(in):: q1(klon, klev) ! specific humidity
36 real, intent(in):: qs1(klon, klev) ! saturation specific humidity
37
38 real, intent(in):: u1(klon, klev), v1(klon, klev)
39 ! zonal wind and meridional velocity (m/s)
40
41 real, intent(in):: p1(klon, klev) ! full level pressure (hPa)
42
43 real, intent(in):: ph1(klon, klev + 1)
44 ! Half level pressure (hPa). These pressures are defined at levels
45 ! intermediate between those of P1, T1, Q1 and QS1. The first
46 ! value of PH should be greater than (i.e. at a lower level than)
47 ! the first value of the array P1.
48
49 integer, intent(out):: iflag1(:) ! (klon)
50 ! Flag for Emanuel conditions.
51
52 ! 0: Moist convection occurs.
53
54 ! 1: Moist convection occurs, but a CFL condition on the
55 ! subsidence warming is violated. This does not cause the scheme
56 ! to terminate.
57
58 ! 2: Moist convection, but no precipitation because ep(inb) < 1e-4
59
60 ! 3: No moist convection because new cbmf is 0 and old cbmf is 0.
61
62 ! 4: No moist convection; atmosphere is not unstable.
63
64 ! 6: No moist convection because ihmin <= minorig.
65
66 ! 7: No moist convection because unreasonable parcel level
67 ! temperature or specific humidity.
68
69 ! 8: No moist convection: lifted condensation level is above the
70 ! 200 mbar level.
71
72 ! 9: No moist convection: cloud base is higher than the level NL-1.
73
74 real, intent(out):: ft1(klon, klev) ! temperature tendency (K/s)
75 real, intent(out):: fq1(klon, klev) ! specific humidity tendency (s-1)
76
77 real, intent(out):: fu1(klon, klev), fv1(klon, klev)
78 ! forcing (tendency) of zonal and meridional velocity (m/s^2)
79
80 real, intent(out):: precip1(klon) ! convective precipitation rate (mm/day)
81
82 real, intent(out):: VPrecip1(klon, klev + 1)
83 ! vertical profile of convective precipitation (kg/m2/s)
84
85 real, intent(inout):: sig1(klon, klev) ! section of adiabatic updraft
86
87 real, intent(inout):: w01(klon, klev)
88 ! vertical velocity within adiabatic updraft
89
90 integer, intent(out):: icb1(klon)
91 integer, intent(inout):: inb1(klon)
92 real, intent(out):: Ma1(klon, klev) ! mass flux of adiabatic updraft
93
94 real, intent(out):: upwd1(klon, klev)
95 ! total upward mass flux (adiabatic + mixed)
96
97 real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed)
98 real, intent(out):: dnwd01(klon, klev) ! unsaturated downward mass flux
99
100 real, intent(out):: qcondc1(klon, klev)
101 ! in-cloud mixing ratio of condensed water
102
103 real, intent(out):: cape1(klon)
104 real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)
105 real, intent(inout):: mp1(klon, klev)
106
107 ! Local:
108
109 real da(klon, klev), phi(klon, klev, klev), mp(klon, klev)
110 integer i, k, il
111 integer nk1(klon)
112 integer icbs1(klon)
113 real plcl1(klon)
114 real tnk1(klon)
115 real qnk1(klon)
116 real gznk1(klon)
117 real pbase1(klon)
118 real buoybase1(klon)
119 real lv1(klon, klev)
120 real cpn1(klon, klev)
121 real tv1(klon, klev)
122 real gz1(klon, klev)
123 real hm1(klon, klev)
124 real h1(klon, klev)
125 real tp1(klon, klev)
126 real tvp1(klon, klev)
127 real clw1(klon, klev)
128 real th1(klon, klev)
129 integer ncum
130
131 ! Compressed fields:
132 integer, allocatable:: idcum(:), iflag(:) ! (ncum)
133 integer nk(klon)
134 integer, allocatable:: icb(:) ! (ncum)
135 integer nent(klon, klev)
136 integer icbs(klon)
137 integer inb(klon)
138 real, allocatable:: plcl(:) ! (ncum)
139 real tnk(klon), qnk(klon), gznk(klon)
140 real t(klon, klev), q(klon, klev), qs(klon, klev)
141 real u(klon, klev), v(klon, klev)
142 real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev)
143 real p(klon, klev) ! pressure at full level, in hPa
144 real ph(klon, klev + 1), tv(klon, klev), tp(klon, klev)
145 real clw(klon, klev)
146 real pbase(klon), buoybase(klon), th(klon, klev)
147 real tvp(klon, klev)
148 real sig(klon, klev), w0(klon, klev)
149 real hp(klon, klev), ep(klon, klev)
150 real buoy(klon, klev)
151 real cape(klon)
152 real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)
153 real uent(klon, klev, klev), vent(klon, klev, klev)
154 real ments(klon, klev, klev), qents(klon, klev, klev)
155 real sij(klon, klev, klev), elij(klon, klev, klev)
156 real qp(klon, klev), up(klon, klev), vp(klon, klev)
157 real wt(klon, klev), water(klon, klev)
158 real, allocatable:: evap(:, :) ! (ncum, nl)
159 real, allocatable:: b(:, :) ! (ncum, nl - 1)
160 real ft(klon, klev), fq(klon, klev)
161 real fu(klon, klev), fv(klon, klev)
162 real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)
163 real Ma(klon, klev), mike(klon, klev), tls(klon, klev)
164 real tps(klon, klev)
165 real precip(klon)
166 real VPrecip(klon, klev + 1)
167 real qcondc(klon, klev) ! cld
168
169 !-------------------------------------------------------------------
170
171 ! SET CONSTANTS AND PARAMETERS
172 CALL cv_thermo
173 CALL cv30_param
174
175 ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS
176
177 do k = 1, klev
178 do i = 1, klon
179 ft1(i, k) = 0.
180 fq1(i, k) = 0.
181 fu1(i, k) = 0.
182 fv1(i, k) = 0.
183 tvp1(i, k) = 0.
184 tp1(i, k) = 0.
185 clw1(i, k) = 0.
186 clw(i, k) = 0.
187 gz1(i, k) = 0.
188 VPrecip1(i, k) = 0.
189 Ma1(i, k) = 0.
190 upwd1(i, k) = 0.
191 dnwd1(i, k) = 0.
192 dnwd01(i, k) = 0.
193 qcondc1(i, k) = 0.
194 end do
195 end do
196
197 precip1 = 0.
198 cape1 = 0.
199 VPrecip1(:, klev + 1) = 0.
200
201 do il = 1, klon
202 sig1(il, klev) = sig1(il, klev) + 1.
203 sig1(il, klev) = min(sig1(il, klev), 12.1)
204 enddo
205
206 CALL cv30_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, &
207 gz1, h1, hm1, th1)
208 CALL cv30_feed(t1, q1, qs1, p1, ph1, gz1, nk1, icb1, iflag1, tnk1, qnk1, &
209 gznk1, plcl1)
210 CALL cv30_undilute1(t1, q1, qs1, gz1, plcl1, p1, nk1, icb1, tp1, tvp1, &
211 clw1, icbs1)
212 CALL cv30_trigger(icb1, plcl1, p1, th1, tv1, tvp1, pbase1, buoybase1, &
213 iflag1, sig1, w01)
214
215 ncum = count(iflag1 == 0)
216
217 IF (ncum > 0) THEN
218 ! Moist convective adjustment is necessary
219 allocate(idcum(ncum), plcl(ncum))
220 allocate(b(ncum, nl - 1), evap(ncum, nl), icb(ncum), iflag(ncum))
221 idcum = pack((/(i, i = 1, klon)/), iflag1 == 0)
222 CALL cv30_compress(iflag1, nk1, icb1, icbs1, plcl1, tnk1, qnk1, gznk1, &
223 pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, cpn1, &
224 p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, nk, icb, icbs, plcl, &
225 tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, th, h, lv, &
226 cpn, p, ph, tv, tp, tvp, clw, sig, w0)
227 CALL cv30_undilute2(icb, icbs(:ncum), nk, tnk, qnk, gznk, t, qs, gz, &
228 p, h, tv, lv, pbase(:ncum), buoybase(:ncum), plcl, inb(:ncum), &
229 tp, tvp, clw, hp, ep, buoy)
230 CALL cv30_closure(icb, inb(:ncum), pbase, p, ph(:ncum, :), tv, buoy, &
231 sig, w0, cape, m)
232 CALL cv30_mixing(icb, nk(:ncum), inb(:ncum), t, q, qs, u, v, h, lv, &
233 hp, ep, clw, m, sig, ment, qent, uent, vent, nent, sij, elij, &
234 ments, qents)
235 CALL cv30_unsat(icb, inb(:ncum), t(:ncum, :nl), q(:ncum, :nl), &
236 qs(:ncum, :nl), gz, u, v, p, ph(:ncum, :), th(:ncum, :nl - 1), &
237 tv, lv, cpn, ep(:ncum, :), clw(:ncum, :), m(:ncum, :), &
238 ment(:ncum, :, :), elij(:ncum, :, :), dtphys, plcl, mp, &
239 qp(:ncum, :nl), up(:ncum, :nl), vp(:ncum, :nl), wt(:ncum, :nl), &
240 water(:ncum, :nl), evap, b)
241 CALL cv30_yield(icb, inb(:ncum), dtphys, t, q, u, v, gz, p, ph, h, hp, &
242 lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp(:ncum, 2:nl), &
243 wt(:ncum, :nl - 1), water(:ncum, :nl), evap, b, ment, qent, uent, &
244 vent, nent, elij, sig, tv, tvp, iflag, precip, VPrecip, ft, fq, &
245 fu, fv, upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc)
246 CALL cv30_tracer(klon, ncum, klev, ment, sij, da, phi)
247 CALL cv30_uncompress(idcum, iflag, precip, VPrecip, sig, w0, ft, fq, &
248 fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, cape, da, phi, mp, &
249 iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, fu1, fv1, inb1, &
250 Ma1, upwd1, dnwd1, dnwd01, qcondc1, cape1, da1, phi1, mp1)
251 ENDIF
252
253 end SUBROUTINE cv_driver
254
255 end module cv_driver_m

  ViewVC Help
Powered by ViewVC 1.1.21