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, & |
8 |
fq1, fu1, fv1, precip1, VPrecip1, cbmf1, sig1, w01, icb1, inb1, delt, & |
9 |
Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, 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 |
! Several versions of convect may be used: |
18 |
! - iflag_con = 3: version lmd |
19 |
! - iflag_con = 4: version 4.3b |
20 |
|
21 |
use clesphys2, only: iflag_con |
22 |
use cv3_compress_m, only: cv3_compress |
23 |
use cv3_feed_m, only: cv3_feed |
24 |
use cv3_mixing_m, only: cv3_mixing |
25 |
use cv3_param_m, only: cv3_param |
26 |
use cv3_prelim_m, only: cv3_prelim |
27 |
use cv3_tracer_m, only: cv3_tracer |
28 |
use cv3_uncompress_m, only: cv3_uncompress |
29 |
use cv3_unsat_m, only: cv3_unsat |
30 |
use cv3_yield_m, only: cv3_yield |
31 |
use cv_feed_m, only: cv_feed |
32 |
use cv_uncompress_m, only: cv_uncompress |
33 |
USE dimphy, ONLY: klev, klon |
34 |
|
35 |
real, intent(in):: t1(klon, klev) ! temperature |
36 |
real, intent(in):: q1(klon, klev) ! specific hum |
37 |
real, intent(in):: qs1(klon, klev) ! sat specific hum |
38 |
real, intent(in):: u1(klon, klev) ! u-wind |
39 |
real, intent(in):: v1(klon, klev) ! v-wind |
40 |
real, intent(in):: p1(klon, klev) ! full level pressure |
41 |
real, intent(in):: ph1(klon, klev + 1) ! half level pressure |
42 |
integer, intent(out):: iflag1(klon) ! flag for Emanuel conditions |
43 |
real, intent(out):: ft1(klon, klev) ! temp tend |
44 |
real, intent(out):: fq1(klon, klev) ! spec hum tend |
45 |
real, intent(out):: fu1(klon, klev) ! u-wind tend |
46 |
real, intent(out):: fv1(klon, klev) ! v-wind tend |
47 |
real, intent(out):: precip1(klon) ! precipitation |
48 |
|
49 |
real, intent(out):: VPrecip1(klon, klev + 1) |
50 |
! vertical profile of precipitation |
51 |
|
52 |
real, intent(inout):: cbmf1(klon) ! cloud base mass flux |
53 |
real, intent(inout):: sig1(klon, klev) ! section adiabatic updraft |
54 |
|
55 |
real, intent(inout):: w01(klon, klev) |
56 |
! vertical velocity within adiabatic updraft |
57 |
|
58 |
integer, intent(out):: icb1(klon) |
59 |
integer, intent(inout):: inb1(klon) |
60 |
real, intent(in):: delt ! time step |
61 |
real Ma1(klon, klev) |
62 |
! Ma1 Real Output mass flux adiabatic updraft |
63 |
|
64 |
real, intent(out):: upwd1(klon, klev) |
65 |
! total upward mass flux (adiab + mixed) |
66 |
|
67 |
real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed) |
68 |
real, intent(out):: dnwd01(klon, klev) ! unsaturated downward mass flux |
69 |
|
70 |
real qcondc1(klon, klev) ! cld |
71 |
! qcondc1 Real Output in-cld mixing ratio of condensed water |
72 |
real wd1(klon) ! gust |
73 |
! wd1 Real Output downdraft velocity scale for sfc fluxes |
74 |
real cape1(klon) |
75 |
! cape1 Real Output CAPE |
76 |
|
77 |
real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev) |
78 |
real, intent(inout):: mp1(klon, klev) |
79 |
|
80 |
! ARGUMENTS |
81 |
|
82 |
! On input: |
83 |
|
84 |
! t: Array of absolute temperature (K) of dimension KLEV, with first |
85 |
! index corresponding to lowest model level. Note that this array |
86 |
! will be altered by the subroutine if dry convective adjustment |
87 |
! occurs and if IPBL is not equal to 0. |
88 |
|
89 |
! q: Array of specific humidity (gm/gm) of dimension KLEV, with first |
90 |
! index corresponding to lowest model level. Must be defined |
91 |
! at same grid levels as T. Note that this array will be altered |
92 |
! if dry convective adjustment occurs and if IPBL is not equal to 0. |
93 |
|
94 |
! qs: Array of saturation specific humidity of dimension KLEV, with first |
95 |
! index corresponding to lowest model level. Must be defined |
96 |
! at same grid levels as T. Note that this array will be altered |
97 |
! if dry convective adjustment occurs and if IPBL is not equal to 0. |
98 |
|
99 |
! u: Array of zonal wind velocity (m/s) of dimension KLEV, witth first |
100 |
! index corresponding with the lowest model level. Defined at |
101 |
! same levels as T. Note that this array will be altered if |
102 |
! dry convective adjustment occurs and if IPBL is not equal to 0. |
103 |
|
104 |
! v: Same as u but for meridional velocity. |
105 |
|
106 |
! p: Array of pressure (mb) of dimension KLEV, with first |
107 |
! index corresponding to lowest model level. Must be defined |
108 |
! at same grid levels as T. |
109 |
|
110 |
! ph: Array of pressure (mb) of dimension KLEV + 1, with first index |
111 |
! corresponding to lowest level. These pressures are defined at |
112 |
! levels intermediate between those of P, T, Q and QS. The first |
113 |
! value of PH should be greater than (i.e. at a lower level than) |
114 |
! the first value of the array P. |
115 |
|
116 |
! nl: The maximum number of levels to which convection can penetrate, plus 1 |
117 |
! NL MUST be less than or equal to KLEV-1. |
118 |
|
119 |
! delt: The model time step (sec) between calls to CONVECT |
120 |
|
121 |
! On Output: |
122 |
|
123 |
! iflag: An output integer whose value denotes the following: |
124 |
! VALUE INTERPRETATION |
125 |
! ----- -------------- |
126 |
! 0 Moist convection occurs. |
127 |
! 1 Moist convection occurs, but a CFL condition |
128 |
! on the subsidence warming is violated. This |
129 |
! does not cause the scheme to terminate. |
130 |
! 2 Moist convection, but no precip because ep(inb) lt 0.0001 |
131 |
! 3 No moist convection because new cbmf is 0 and old cbmf is 0. |
132 |
! 4 No moist convection; atmosphere is not |
133 |
! unstable |
134 |
! 6 No moist convection because ihmin le minorig. |
135 |
! 7 No moist convection because unreasonable |
136 |
! parcel level temperature or specific humidity. |
137 |
! 8 No moist convection: lifted condensation |
138 |
! level is above the 200 mb level. |
139 |
! 9 No moist convection: cloud base is higher |
140 |
! then the level NL-1. |
141 |
|
142 |
! ft: Array of temperature tendency (K/s) of dimension KLEV, defined at same |
143 |
! grid levels as T, Q, QS and P. |
144 |
|
145 |
! fq: Array of specific humidity tendencies ((gm/gm)/s) of dimension KLEV, |
146 |
! defined at same grid levels as T, Q, QS and P. |
147 |
|
148 |
! fu: Array of forcing of zonal velocity (m/s^2) of dimension KLEV, |
149 |
! defined at same grid levels as T. |
150 |
|
151 |
! fv: Same as FU, but for forcing of meridional velocity. |
152 |
|
153 |
! precip: Scalar convective precipitation rate (mm/day). |
154 |
|
155 |
! VPrecip: Vertical profile of convective precipitation (kg/m2/s). |
156 |
|
157 |
! wd: A convective downdraft velocity scale. For use in surface |
158 |
! flux parameterizations. See convect.ps file for details. |
159 |
|
160 |
! tprime: A convective downdraft temperature perturbation scale (K). |
161 |
! For use in surface flux parameterizations. See convect.ps |
162 |
! file for details. |
163 |
|
164 |
! qprime: A convective downdraft specific humidity |
165 |
! perturbation scale (gm/gm). |
166 |
! For use in surface flux parameterizations. See convect.ps |
167 |
! file for details. |
168 |
|
169 |
! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST |
170 |
! BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT |
171 |
! ITS NEXT CALL. That is, the value of CBMF must be "remembered" |
172 |
! by the calling program between calls to CONVECT. |
173 |
|
174 |
! det: Array of detrainment mass flux of dimension KLEV. |
175 |
|
176 |
! Local arrays |
177 |
|
178 |
real da(klon, klev), phi(klon, klev, klev), mp(klon, klev) |
179 |
|
180 |
integer i, k, il |
181 |
integer icbmax |
182 |
integer nk1(klon) |
183 |
integer icbs1(klon) |
184 |
|
185 |
real plcl1(klon) |
186 |
real tnk1(klon) |
187 |
real qnk1(klon) |
188 |
real gznk1(klon) |
189 |
real pbase1(klon) |
190 |
real buoybase1(klon) |
191 |
|
192 |
real lv1(klon, klev) |
193 |
real cpn1(klon, klev) |
194 |
real tv1(klon, klev) |
195 |
real gz1(klon, klev) |
196 |
real hm1(klon, klev) |
197 |
real h1(klon, klev) |
198 |
real tp1(klon, klev) |
199 |
real tvp1(klon, klev) |
200 |
real clw1(klon, klev) |
201 |
real th1(klon, klev) |
202 |
|
203 |
integer ncum |
204 |
|
205 |
! (local) compressed fields: |
206 |
|
207 |
integer idcum(klon) |
208 |
integer iflag(klon), nk(klon), icb(klon) |
209 |
integer nent(klon, klev) |
210 |
integer icbs(klon) |
211 |
integer inb(klon), inbis(klon) |
212 |
|
213 |
real cbmf(klon), plcl(klon), tnk(klon), qnk(klon), gznk(klon) |
214 |
real t(klon, klev), q(klon, klev), qs(klon, klev) |
215 |
real u(klon, klev), v(klon, klev) |
216 |
real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev) |
217 |
real p(klon, klev), ph(klon, klev + 1), tv(klon, klev), tp(klon, klev) |
218 |
real clw(klon, klev) |
219 |
real dph(klon, klev) |
220 |
real pbase(klon), buoybase(klon), th(klon, klev) |
221 |
real tvp(klon, klev) |
222 |
real sig(klon, klev), w0(klon, klev) |
223 |
real hp(klon, klev), ep(klon, klev), sigp(klon, klev) |
224 |
real frac(klon), buoy(klon, klev) |
225 |
real cape(klon) |
226 |
real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev) |
227 |
real uent(klon, klev, klev), vent(klon, klev, klev) |
228 |
real ments(klon, klev, klev), qents(klon, klev, klev) |
229 |
real sij(klon, klev, klev), elij(klon, klev, klev) |
230 |
real qp(klon, klev), up(klon, klev), vp(klon, klev) |
231 |
real wt(klon, klev), water(klon, klev), evap(klon, klev) |
232 |
real b(klon, klev), ft(klon, klev), fq(klon, klev) |
233 |
real fu(klon, klev), fv(klon, klev) |
234 |
real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev) |
235 |
real Ma(klon, klev), mike(klon, klev), tls(klon, klev) |
236 |
real tps(klon, klev), qprime(klon), tprime(klon) |
237 |
real precip(klon) |
238 |
real VPrecip(klon, klev + 1) |
239 |
real qcondc(klon, klev) ! cld |
240 |
real wd(klon) ! gust |
241 |
|
242 |
!------------------------------------------------------------------- |
243 |
|
244 |
! SET CONSTANTS AND PARAMETERS |
245 |
|
246 |
! set simulation flags: |
247 |
! (common cvflag) |
248 |
|
249 |
CALL cv_flag |
250 |
|
251 |
! set thermodynamical constants: |
252 |
! (common cvthermo) |
253 |
|
254 |
CALL cv_thermo |
255 |
|
256 |
! set convect parameters |
257 |
|
258 |
! includes microphysical parameters and parameters that |
259 |
! control the rate of approach to quasi-equilibrium) |
260 |
! (common cvparam) |
261 |
|
262 |
if (iflag_con == 3) CALL cv3_param(klev, delt) |
263 |
|
264 |
! INITIALIZE OUTPUT ARRAYS AND PARAMETERS |
265 |
|
266 |
do k = 1, klev |
267 |
do i = 1, klon |
268 |
ft1(i, k) = 0.0 |
269 |
fq1(i, k) = 0.0 |
270 |
fu1(i, k) = 0.0 |
271 |
fv1(i, k) = 0.0 |
272 |
tvp1(i, k) = 0.0 |
273 |
tp1(i, k) = 0.0 |
274 |
clw1(i, k) = 0.0 |
275 |
!ym |
276 |
clw(i, k) = 0.0 |
277 |
gz1(i, k) = 0. |
278 |
VPrecip1(i, k) = 0. |
279 |
Ma1(i, k) = 0.0 |
280 |
upwd1(i, k) = 0.0 |
281 |
dnwd1(i, k) = 0.0 |
282 |
dnwd01(i, k) = 0.0 |
283 |
qcondc1(i, k) = 0.0 |
284 |
end do |
285 |
end do |
286 |
|
287 |
do i = 1, klon |
288 |
precip1(i) = 0.0 |
289 |
iflag1(i) = 0 |
290 |
wd1(i) = 0.0 |
291 |
cape1(i) = 0.0 |
292 |
VPrecip1(i, klev + 1) = 0.0 |
293 |
end do |
294 |
|
295 |
if (iflag_con == 3) then |
296 |
do il = 1, klon |
297 |
sig1(il, klev) = sig1(il, klev) + 1. |
298 |
sig1(il, klev) = min(sig1(il, klev), 12.1) |
299 |
enddo |
300 |
endif |
301 |
|
302 |
! CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY |
303 |
|
304 |
if (iflag_con == 3) then |
305 |
CALL cv3_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, & |
306 |
gz1, h1, hm1, th1) |
307 |
else |
308 |
! iflag_con == 4 |
309 |
CALL cv_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, & |
310 |
gz1, h1, hm1) |
311 |
endif |
312 |
|
313 |
! CONVECTIVE FEED |
314 |
|
315 |
if (iflag_con == 3) then |
316 |
CALL cv3_feed(klon, klev, t1, q1, qs1, p1, ph1, gz1, nk1, icb1, & |
317 |
icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! klev->na |
318 |
else |
319 |
! iflag_con == 4 |
320 |
CALL cv_feed(klon, klev, t1, q1, qs1, p1, hm1, gz1, nk1, icb1, icbmax, & |
321 |
iflag1, tnk1, qnk1, gznk1, plcl1) |
322 |
endif |
323 |
|
324 |
! UNDILUTE (ADIABATIC) UPDRAFT / 1st part |
325 |
! (up through ICB for convect4, up through ICB + 1 for convect3) |
326 |
! Calculates the lifted parcel virtual temperature at nk, the |
327 |
! actual temperature, and the adiabatic liquid water content. |
328 |
|
329 |
if (iflag_con == 3) then |
330 |
CALL cv3_undilute1(klon, klev, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1, & |
331 |
tp1, tvp1, clw1, icbs1) ! klev->na |
332 |
else |
333 |
! iflag_con == 4 |
334 |
CALL cv_undilute1(klon, klev, t1, q1, qs1, gz1, p1, nk1, icb1, icbmax, & |
335 |
tp1, tvp1, clw1) |
336 |
endif |
337 |
|
338 |
! TRIGGERING |
339 |
|
340 |
if (iflag_con == 3) then |
341 |
CALL cv3_trigger(klon, klev, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, & |
342 |
buoybase1, iflag1, sig1, w01) ! klev->na |
343 |
else |
344 |
! iflag_con == 4 |
345 |
CALL cv_trigger(klon, klev, icb1, cbmf1, tv1, tvp1, iflag1) |
346 |
end if |
347 |
|
348 |
! Moist convective adjustment is necessary |
349 |
|
350 |
ncum = 0 |
351 |
do i = 1, klon |
352 |
if (iflag1(i) == 0) then |
353 |
ncum = ncum + 1 |
354 |
idcum(ncum) = i |
355 |
endif |
356 |
end do |
357 |
|
358 |
IF (ncum > 0) THEN |
359 |
! COMPRESS THE FIELDS |
360 |
! (-> vectorization over convective gridpoints) |
361 |
|
362 |
if (iflag_con == 3) then |
363 |
CALL cv3_compress(klon, klon, ncum, klev, iflag1, nk1, icb1, icbs1, & |
364 |
plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, qs1, u1, & |
365 |
v1, gz1, th1, h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1, & |
366 |
sig1, w01, iflag, nk, icb, icbs, plcl, tnk, qnk, gznk, pbase, & |
367 |
buoybase, t, q, qs, u, v, gz, th, h, lv, cpn, p, ph, tv, tp, & |
368 |
tvp, clw, sig, w0) |
369 |
else |
370 |
! iflag_con == 4 |
371 |
CALL cv_compress(klon, klon, ncum, klev, iflag1, nk1, icb1, cbmf1, & |
372 |
plcl1, tnk1, qnk1, gznk1, t1, q1, qs1, u1, v1, gz1, h1, lv1, & |
373 |
cpn1, p1, ph1, tv1, tp1, tvp1, clw1, iflag, nk, icb, cbmf, & |
374 |
plcl, tnk, qnk, gznk, t, q, qs, u, v, gz, h, lv, cpn, p, ph, & |
375 |
tv, tp, tvp, clw, dph) |
376 |
endif |
377 |
|
378 |
! UNDILUTE (ADIABATIC) UPDRAFT / second part : |
379 |
! FIND THE REST OF THE LIFTED PARCEL TEMPERATURES |
380 |
! & |
381 |
! COMPUTE THE PRECIPITATION EFFICIENCIES AND THE |
382 |
! FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD |
383 |
! & |
384 |
! FIND THE LEVEL OF NEUTRAL BUOYANCY |
385 |
|
386 |
if (iflag_con == 3) then |
387 |
CALL cv3_undilute2(klon, ncum, klev, icb, icbs, nk, tnk, qnk, gznk, & |
388 |
t, qs, gz, p, h, tv, lv, pbase, buoybase, plcl, inb, tp, & |
389 |
tvp, clw, hp, ep, sigp, buoy) !na->klev |
390 |
else |
391 |
! iflag_con == 4 |
392 |
CALL cv_undilute2(klon, ncum, klev, icb, nk, tnk, qnk, gznk, t, & |
393 |
qs, gz, p, dph, h, tv, lv, inb, inbis, tp, tvp, clw, hp, ep, & |
394 |
sigp, frac) |
395 |
endif |
396 |
|
397 |
! CLOSURE |
398 |
|
399 |
if (iflag_con == 3) then |
400 |
CALL cv3_closure(klon, ncum, klev, icb, inb, pbase, p, ph, tv, & |
401 |
buoy, sig, w0, cape, m) ! na->klev |
402 |
else |
403 |
! iflag_con == 4 |
404 |
CALL cv_closure(klon, ncum, klev, nk, icb, tv, tvp, p, ph, dph, & |
405 |
plcl, cpn, iflag, cbmf) |
406 |
endif |
407 |
|
408 |
! MIXING |
409 |
|
410 |
if (iflag_con == 3) then |
411 |
CALL cv3_mixing(klon, ncum, klev, klev, icb, nk, inb, t, q, qs, u, & |
412 |
v, h, lv, hp, ep, clw, m, sig, ment, qent, uent, vent, nent, & |
413 |
sij, elij, ments, qents) |
414 |
else |
415 |
! iflag_con == 4 |
416 |
CALL cv_mixing(klon, ncum, klev, icb, nk, inb, inbis, ph, t, q, qs, & |
417 |
u, v, h, lv, qnk, hp, tv, tvp, ep, clw, cbmf, m, ment, qent, & |
418 |
uent, vent, nent, sij, elij) |
419 |
endif |
420 |
|
421 |
! UNSATURATED (PRECIPITATING) DOWNDRAFTS |
422 |
|
423 |
if (iflag_con == 3) then |
424 |
CALL cv3_unsat(klon, ncum, klev, klev, icb, inb, t, q, qs, gz, u, & |
425 |
v, p, ph, th, tv, lv, cpn, ep, sigp, clw, m, ment, elij, delt, & |
426 |
plcl, mp, qp, up, vp, wt, water, evap, b)! na->klev |
427 |
else |
428 |
! iflag_con == 4 |
429 |
CALL cv_unsat(klon, ncum, klev, inb, t, q, qs, gz, u, v, p, ph, h, & |
430 |
lv, ep, sigp, clw, m, ment, elij, iflag, mp, qp, up, vp, wt, & |
431 |
water, evap) |
432 |
endif |
433 |
|
434 |
! YIELD |
435 |
! (tendencies, precipitation, variables of interface with other |
436 |
! processes, etc) |
437 |
|
438 |
if (iflag_con == 3) then |
439 |
CALL cv3_yield(klon, ncum, klev, klev, icb, inb, delt, t, q, u, v, & |
440 |
gz, p, ph, h, hp, lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp, & |
441 |
wt, water, evap, b, ment, qent, uent, vent, nent, elij, sig, & |
442 |
tv, tvp, iflag, precip, VPrecip, ft, fq, fu, fv, upwd, dnwd, & |
443 |
dnwd0, ma, mike, tls, tps, qcondc, wd)! na->klev |
444 |
else |
445 |
! iflag_con == 4 |
446 |
CALL cv_yield(klon, ncum, klev, nk, icb, inb, delt, t, q, u, v, gz, & |
447 |
p, ph, h, hp, lv, cpn, ep, clw, frac, m, mp, qp, up, vp, wt, & |
448 |
water, evap, ment, qent, uent, vent, nent, elij, tv, tvp, & |
449 |
iflag, wd, qprime, tprime, precip, cbmf, ft, fq, fu, fv, Ma, & |
450 |
qcondc) |
451 |
endif |
452 |
|
453 |
! passive tracers |
454 |
|
455 |
if (iflag_con == 3) CALL cv3_tracer(klon, ncum, klev, ment, sij, da, phi) |
456 |
|
457 |
! UNCOMPRESS THE FIELDS |
458 |
|
459 |
! set iflag1 = 42 for non convective points |
460 |
do i = 1, klon |
461 |
iflag1(i) = 42 |
462 |
end do |
463 |
|
464 |
if (iflag_con == 3) then |
465 |
CALL cv3_uncompress(idcum(:ncum), iflag, precip, VPrecip, sig, w0, & |
466 |
ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, wd, cape, & |
467 |
da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, & |
468 |
fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, & |
469 |
cape1, da1, phi1, mp1) |
470 |
else |
471 |
! iflag_con == 4 |
472 |
CALL cv_uncompress(idcum(:ncum), iflag, precip, cbmf, ft, fq, fu, & |
473 |
fv, Ma, qcondc, iflag1, precip1, cbmf1, ft1, fq1, fu1, fv1, & |
474 |
Ma1, qcondc1) |
475 |
endif |
476 |
ENDIF ! ncum>0 |
477 |
|
478 |
end SUBROUTINE cv_driver |
479 |
|
480 |
end module cv_driver_m |