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, delt, 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 cv30_closure_m, only: cv30_closure |
18 |
use cv30_compress_m, only: cv30_compress |
19 |
use cv30_feed_m, only: cv30_feed |
20 |
use cv30_mixing_m, only: cv30_mixing |
21 |
use cv30_param_m, only: cv30_param, nl |
22 |
use cv30_prelim_m, only: cv30_prelim |
23 |
use cv30_tracer_m, only: cv30_tracer |
24 |
use cv30_trigger_m, only: cv30_trigger |
25 |
use cv30_uncompress_m, only: cv30_uncompress |
26 |
use cv30_undilute2_m, only: cv30_undilute2 |
27 |
use cv30_unsat_m, only: cv30_unsat |
28 |
use cv30_yield_m, only: cv30_yield |
29 |
use cv_thermo_m, only: cv_thermo |
30 |
USE dimphy, ONLY: klev, klon |
31 |
|
32 |
real, intent(in):: t1(klon, klev) ! temperature (K) |
33 |
real, intent(in):: q1(klon, klev) ! specific humidity |
34 |
real, intent(in):: qs1(klon, klev) ! saturation specific humidity |
35 |
|
36 |
real, intent(in):: u1(klon, klev), v1(klon, klev) |
37 |
! zonal wind and meridional velocity (m/s) |
38 |
|
39 |
real, intent(in):: p1(klon, klev) ! full level pressure (hPa) |
40 |
|
41 |
real, intent(in):: ph1(klon, klev + 1) |
42 |
! Half level pressure (hPa). These pressures are defined at levels |
43 |
! intermediate between those of P1, T1, Q1 and QS1. The first |
44 |
! value of PH should be greater than (i.e. at a lower level than) |
45 |
! the first value of the array P1. |
46 |
|
47 |
integer, intent(out):: iflag1(klon) |
48 |
! Flag for Emanuel conditions. |
49 |
|
50 |
! 0: Moist convection occurs. |
51 |
|
52 |
! 1: Moist convection occurs, but a CFL condition on the |
53 |
! subsidence warming is violated. This does not cause the scheme |
54 |
! to terminate. |
55 |
|
56 |
! 2: Moist convection, but no precipitation because ep(inb) < 1e-4 |
57 |
|
58 |
! 3: No moist convection because new cbmf is 0 and old cbmf is 0. |
59 |
|
60 |
! 4: No moist convection; atmosphere is not unstable |
61 |
|
62 |
! 6: No moist convection because ihmin le minorig. |
63 |
|
64 |
! 7: No moist convection because unreasonable parcel level |
65 |
! temperature or specific humidity. |
66 |
|
67 |
! 8: No moist convection: lifted condensation level is above the |
68 |
! 200 mb level. |
69 |
|
70 |
! 9: No moist convection: cloud base is higher then the level NL-1. |
71 |
|
72 |
real, intent(out):: ft1(klon, klev) ! temperature tendency (K/s) |
73 |
real, intent(out):: fq1(klon, klev) ! specific humidity tendency (s-1) |
74 |
|
75 |
real, intent(out):: fu1(klon, klev), fv1(klon, klev) |
76 |
! forcing (tendency) of zonal and meridional velocity (m/s^2) |
77 |
|
78 |
real, intent(out):: precip1(klon) ! convective precipitation rate (mm/day) |
79 |
|
80 |
real, intent(out):: VPrecip1(klon, klev + 1) |
81 |
! vertical profile of convective precipitation (kg/m2/s) |
82 |
|
83 |
real, intent(inout):: sig1(klon, klev) ! section of adiabatic updraft |
84 |
|
85 |
real, intent(inout):: w01(klon, klev) |
86 |
! vertical velocity within adiabatic updraft |
87 |
|
88 |
integer, intent(out):: icb1(klon) |
89 |
integer, intent(inout):: inb1(klon) |
90 |
real, intent(in):: delt ! the model time step (sec) between calls |
91 |
|
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 icbmax |
112 |
integer nk1(klon) |
113 |
integer icbs1(klon) |
114 |
real plcl1(klon) |
115 |
real tnk1(klon) |
116 |
real qnk1(klon) |
117 |
real gznk1(klon) |
118 |
real pbase1(klon) |
119 |
real buoybase1(klon) |
120 |
real lv1(klon, klev) |
121 |
real cpn1(klon, klev) |
122 |
real tv1(klon, klev) |
123 |
real gz1(klon, klev) |
124 |
real hm1(klon, klev) |
125 |
real h1(klon, klev) |
126 |
real tp1(klon, klev) |
127 |
real tvp1(klon, klev) |
128 |
real clw1(klon, klev) |
129 |
real th1(klon, klev) |
130 |
integer ncum |
131 |
|
132 |
! Compressed fields: |
133 |
integer idcum(klon) |
134 |
integer iflag(klon), nk(klon), icb(klon) |
135 |
integer nent(klon, klev) |
136 |
integer icbs(klon) |
137 |
integer inb(klon) |
138 |
real plcl(klon), tnk(klon), qnk(klon), gznk(klon) |
139 |
real t(klon, klev), q(klon, klev), qs(klon, klev) |
140 |
real u(klon, klev), v(klon, klev) |
141 |
real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev) |
142 |
real p(klon, klev), ph(klon, klev + 1), tv(klon, klev), tp(klon, klev) |
143 |
real clw(klon, klev) |
144 |
real pbase(klon), buoybase(klon), th(klon, klev) |
145 |
real tvp(klon, klev) |
146 |
real sig(klon, klev), w0(klon, klev) |
147 |
real hp(klon, klev), ep(klon, klev), sigp(klon, klev) |
148 |
real buoy(klon, klev) |
149 |
real cape(klon) |
150 |
real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev) |
151 |
real uent(klon, klev, klev), vent(klon, klev, klev) |
152 |
real ments(klon, klev, klev), qents(klon, klev, klev) |
153 |
real sij(klon, klev, klev), elij(klon, klev, klev) |
154 |
real qp(klon, klev), up(klon, klev), vp(klon, klev) |
155 |
real wt(klon, klev), water(klon, klev), evap(klon, klev) |
156 |
real, allocatable:: b(:, :) ! (ncum, nl - 1) |
157 |
real ft(klon, klev), fq(klon, klev) |
158 |
real fu(klon, klev), fv(klon, klev) |
159 |
real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev) |
160 |
real Ma(klon, klev), mike(klon, klev), tls(klon, klev) |
161 |
real tps(klon, klev) |
162 |
real precip(klon) |
163 |
real VPrecip(klon, klev + 1) |
164 |
real qcondc(klon, klev) ! cld |
165 |
|
166 |
!------------------------------------------------------------------- |
167 |
|
168 |
! SET CONSTANTS AND PARAMETERS |
169 |
|
170 |
! set thermodynamical constants: |
171 |
CALL cv_thermo |
172 |
|
173 |
! set convect parameters |
174 |
! includes microphysical parameters and parameters that |
175 |
! control the rate of approach to quasi-equilibrium) |
176 |
! (common cvparam) |
177 |
CALL cv30_param(delt) |
178 |
|
179 |
! INITIALIZE OUTPUT ARRAYS AND PARAMETERS |
180 |
|
181 |
do k = 1, klev |
182 |
do i = 1, klon |
183 |
ft1(i, k) = 0. |
184 |
fq1(i, k) = 0. |
185 |
fu1(i, k) = 0. |
186 |
fv1(i, k) = 0. |
187 |
tvp1(i, k) = 0. |
188 |
tp1(i, k) = 0. |
189 |
clw1(i, k) = 0. |
190 |
clw(i, k) = 0. |
191 |
gz1(i, k) = 0. |
192 |
VPrecip1(i, k) = 0. |
193 |
Ma1(i, k) = 0. |
194 |
upwd1(i, k) = 0. |
195 |
dnwd1(i, k) = 0. |
196 |
dnwd01(i, k) = 0. |
197 |
qcondc1(i, k) = 0. |
198 |
end do |
199 |
end do |
200 |
|
201 |
do i = 1, klon |
202 |
precip1(i) = 0. |
203 |
iflag1(i) = 0 |
204 |
cape1(i) = 0. |
205 |
VPrecip1(i, klev + 1) = 0. |
206 |
end do |
207 |
|
208 |
do il = 1, klon |
209 |
sig1(il, klev) = sig1(il, klev) + 1. |
210 |
sig1(il, klev) = min(sig1(il, klev), 12.1) |
211 |
enddo |
212 |
|
213 |
! CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY |
214 |
CALL cv30_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, & |
215 |
gz1, h1, hm1, th1) |
216 |
|
217 |
! CONVECTIVE FEED |
218 |
CALL cv30_feed(klon, klev, t1, q1, qs1, p1, ph1, gz1, nk1, icb1, & |
219 |
icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! klev->na |
220 |
|
221 |
CALL cv30_undilute1(klon, klev, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1, & |
222 |
tp1, tvp1, clw1, icbs1) ! klev->na |
223 |
|
224 |
! TRIGGERING |
225 |
CALL cv30_trigger(klon, klev, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, & |
226 |
buoybase1, iflag1, sig1, w01) ! klev->na |
227 |
|
228 |
! Moist convective adjustment is necessary |
229 |
|
230 |
ncum = 0 |
231 |
do i = 1, klon |
232 |
if (iflag1(i) == 0) then |
233 |
ncum = ncum + 1 |
234 |
idcum(ncum) = i |
235 |
endif |
236 |
end do |
237 |
|
238 |
IF (ncum > 0) THEN |
239 |
allocate(b(ncum, nl - 1)) |
240 |
CALL cv30_compress(ncum, iflag1, nk1, icb1, icbs1, plcl1, tnk1, qnk1, & |
241 |
gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, & |
242 |
cpn1, p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, iflag, nk, icb, & |
243 |
icbs, plcl, tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, & |
244 |
th, h, lv, cpn, p, ph, tv, tp, tvp, clw, sig, w0) |
245 |
CALL cv30_undilute2(ncum, icb, icbs, nk, tnk, qnk, gznk, t, qs, gz, p, & |
246 |
h, tv, lv, pbase, buoybase, plcl, inb(:ncum), tp, tvp, clw, hp, & |
247 |
ep, sigp, buoy) |
248 |
|
249 |
! CLOSURE |
250 |
CALL cv30_closure(klon, ncum, klev, icb, inb, pbase, p, ph, tv, & |
251 |
buoy, sig, w0, cape, m) ! na->klev |
252 |
|
253 |
! MIXING |
254 |
CALL cv30_mixing(klon, ncum, klev, klev, icb, nk, inb, t, q, qs, u, & |
255 |
v, h, lv, hp, ep, clw, m, sig, ment, qent, uent, vent, nent, & |
256 |
sij, elij, ments, qents) |
257 |
|
258 |
! Unsaturated (precipitating) downdrafts |
259 |
CALL cv30_unsat(icb(:ncum), inb(:ncum), t, q, qs, gz, u, v, p, ph, th, & |
260 |
tv, lv, cpn, ep(:ncum, :), sigp(:ncum, :), clw(:ncum, :), & |
261 |
m(:ncum, :), ment(:ncum, :, :), elij(:ncum, :, :), delt, plcl, & |
262 |
mp, qp(:ncum, :nl), up(:ncum, :nl), vp(:ncum, :nl), wt, water, & |
263 |
evap, b) |
264 |
|
265 |
! Yield (tendencies, precipitation, variables of interface with |
266 |
! other processes, etc) |
267 |
CALL cv30_yield(icb(:ncum), inb(:ncum), delt, t, q, u, v, gz, p, ph, & |
268 |
h, hp, lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp, wt, & |
269 |
water(:ncum, :nl), evap(:ncum, :nl), b, ment, qent, uent, vent, & |
270 |
nent, elij, sig, tv, tvp, iflag, precip, VPrecip, ft, fq, fu, fv, & |
271 |
upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc) |
272 |
|
273 |
CALL cv30_tracer(klon, ncum, klev, ment, sij, da, phi) |
274 |
|
275 |
! UNCOMPRESS THE FIELDS |
276 |
iflag1 = 42 ! for non convective points |
277 |
CALL cv30_uncompress(idcum(:ncum), iflag, precip, VPrecip, sig, w0, & |
278 |
ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, cape, & |
279 |
da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, & |
280 |
fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, cape1, da1, & |
281 |
phi1, mp1) |
282 |
ENDIF |
283 |
|
284 |
end SUBROUTINE cv_driver |
285 |
|
286 |
end module cv_driver_m |