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

Annotation of /trunk/Sources/phylmd/cv_driver.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 91 - (hide annotations)
Wed Mar 26 17:18:58 2014 UTC (10 years, 1 month ago) by guez
Original Path: trunk/phylmd/cv_driver.f
File size: 23398 byte(s)
Removed unused variables lock_startdate and time_stamp of module
calendar.

Noticed that physiq does not change the surface pressure. So removed
arguments ps and dpfi of subroutine addfi. dpfi was always 0. The
computation of ps in addfi included some averaging at the poles. In
principle, this does not change ps but in practice it does because of
finite numerical precision. So the results of the simulation are
changed. Removed arguments ps and dpfi of calfis. Removed argument
d_ps of physiq.

du at the poles is not computed by dudv1, so declare only the
corresponding latitudes in dudv1. caldyn passes only a section of the
array dudyn as argument.

Removed variable niadv of module iniadvtrac_m.

Declared arguments of exner_hyb as assumed-shape arrays and made all
other horizontal sizes in exner_hyb dynamic. This allows the external
program test_disvert to use exner_hyb at a single horizontal position.

1 guez 52 module cv_driver_m
2 guez 3
3 guez 52 implicit none
4 guez 3
5 guez 52 contains
6 guez 3
7 guez 69 SUBROUTINE cv_driver(len, nd, ndp1, ntra, t1, q1, qs1, u1, v1, tra1, p1, &
8     ph1, iflag1, ft1, fq1, fu1, fv1, ftra1, precip1, VPrecip1, cbmf1, &
9     sig1, w01, icb1, inb1, delt, Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, &
10     cape1, da1, phi1, mp1)
11 guez 3
12 guez 91 ! From LMDZ4/libf/phylmd/cv_driver.F, version 1.3, 2005/04/15 12:36:17
13 guez 3
14 guez 69 ! Main driver for convection
15    
16 guez 91 ! S. Bony, March 2002:
17 guez 72
18     ! Several modules corresponding to different physical processes
19    
20     ! Several versions of convect may be used:
21 guez 91 ! - iflag_con = 3: version lmd (previously named convect3)
22     ! - iflag_con = 4: version 4.3b (vect. version, previously convect1/2)
23 guez 72
24     ! Plus tard :
25 guez 91 ! - iflag_con = 5: version lmd with ice (previously named convectg)
26 guez 72
27     ! S. Bony, Oct 2002:
28     ! Vectorization of convect3 (ie version lmd)
29    
30 guez 69 use clesphys2, only: iflag_con
31 guez 91 use cv3_compress_m, only: cv3_compress
32 guez 69 use cv3_param_m, only: cv3_param
33 guez 62 USE dimphy, ONLY: klev, klon
34    
35 guez 52 ! PARAMETERS:
36     ! Name Type Usage Description
37     ! ---------- ---------- ------- ----------------------------
38 guez 62
39 guez 52 ! len Integer Input first (i) dimension
40     ! nd Integer Input vertical (k) dimension
41     ! ndp1 Integer Input nd + 1
42     ! ntra Integer Input number of tracors
43     ! t1 Real Input temperature
44     ! q1 Real Input specific hum
45     ! qs1 Real Input sat specific hum
46     ! u1 Real Input u-wind
47     ! v1 Real Input v-wind
48     ! tra1 Real Input tracors
49     ! p1 Real Input full level pressure
50     ! ph1 Real Input half level pressure
51     ! iflag1 Integer Output flag for Emanuel conditions
52     ! ft1 Real Output temp tend
53     ! fq1 Real Output spec hum tend
54     ! fu1 Real Output u-wind tend
55     ! fv1 Real Output v-wind tend
56     ! ftra1 Real Output tracor tend
57     ! precip1 Real Output precipitation
58     ! VPrecip1 Real Output vertical profile of precipitations
59     ! cbmf1 Real Output cloud base mass flux
60     ! delt Real Input time step
61     ! Ma1 Real Output mass flux adiabatic updraft
62     ! qcondc1 Real Output in-cld mixing ratio of condensed water
63     ! wd1 Real Output downdraft velocity scale for sfc fluxes
64     ! cape1 Real Output CAPE
65 guez 62
66 guez 52 integer len
67     integer nd
68     integer ndp1
69 guez 69 integer, intent(in):: ntra
70 guez 52 real, intent(in):: t1(len, nd)
71     real q1(len, nd)
72     real qs1(len, nd)
73 guez 91 real, intent(in):: u1(len, nd)
74     real, intent(in):: v1(len, nd)
75 guez 72 real, intent(in):: tra1(len, nd, ntra)
76 guez 52 real p1(len, nd)
77     real ph1(len, ndp1)
78     integer iflag1(len)
79     real ft1(len, nd)
80     real fq1(len, nd)
81     real fu1(len, nd)
82     real fv1(len, nd)
83 guez 72 real ftra1(len, nd, ntra)
84 guez 52 real precip1(len)
85 guez 72 real VPrecip1(len, nd+1)
86 guez 52 real cbmf1(len)
87 guez 72 real, intent(inout):: sig1(klon, klev) ! section adiabatic updraft
88    
89     real, intent(inout):: w01(klon, klev)
90     ! vertical velocity within adiabatic updraft
91    
92     integer icb1(klon)
93     integer inb1(klon)
94     real, intent(in):: delt
95 guez 52 real Ma1(len, nd)
96 guez 62 real, intent(out):: upwd1(len, nd) ! total upward mass flux (adiab+mixed)
97     real, intent(out):: dnwd1(len, nd) ! saturated downward mass flux (mixed)
98     real, intent(out):: dnwd01(len, nd) ! unsaturated downward mass flux
99 guez 3
100 guez 52 real qcondc1(len, nd) ! cld
101     real wd1(len) ! gust
102     real cape1(len)
103 guez 3
104 guez 52 real da1(len, nd), phi1(len, nd, nd), mp1(len, nd)
105 guez 3
106 guez 52 !-------------------------------------------------------------------
107     ! --- ARGUMENTS
108     !-------------------------------------------------------------------
109     ! --- On input:
110 guez 62
111 guez 52 ! t: Array of absolute temperature (K) of dimension ND, with first
112     ! index corresponding to lowest model level. Note that this array
113     ! will be altered by the subroutine if dry convective adjustment
114     ! occurs and if IPBL is not equal to 0.
115 guez 62
116 guez 52 ! q: Array of specific humidity (gm/gm) of dimension ND, with first
117     ! index corresponding to lowest model level. Must be defined
118     ! at same grid levels as T. Note that this array will be altered
119     ! if dry convective adjustment occurs and if IPBL is not equal to 0.
120 guez 62
121 guez 52 ! qs: Array of saturation specific humidity of dimension ND, with first
122     ! index corresponding to lowest model level. Must be defined
123     ! at same grid levels as T. Note that this array will be altered
124     ! if dry convective adjustment occurs and if IPBL is not equal to 0.
125 guez 62
126 guez 52 ! u: Array of zonal wind velocity (m/s) of dimension ND, witth first
127     ! index corresponding with the lowest model level. Defined at
128     ! same levels as T. Note that this array will be altered if
129     ! dry convective adjustment occurs and if IPBL is not equal to 0.
130 guez 62
131 guez 52 ! v: Same as u but for meridional velocity.
132 guez 62
133 guez 52 ! tra: Array of passive tracer mixing ratio, of dimensions (ND, NTRA),
134     ! where NTRA is the number of different tracers. If no
135     ! convective tracer transport is needed, define a dummy
136     ! input array of dimension (ND, 1). Tracers are defined at
137     ! same vertical levels as T. Note that this array will be altered
138     ! if dry convective adjustment occurs and if IPBL is not equal to 0.
139 guez 62
140 guez 52 ! p: Array of pressure (mb) of dimension ND, with first
141     ! index corresponding to lowest model level. Must be defined
142     ! at same grid levels as T.
143 guez 62
144 guez 52 ! ph: Array of pressure (mb) of dimension ND+1, with first index
145     ! corresponding to lowest level. These pressures are defined at
146     ! levels intermediate between those of P, T, Q and QS. The first
147     ! value of PH should be greater than (i.e. at a lower level than)
148     ! the first value of the array P.
149 guez 62
150 guez 52 ! nl: The maximum number of levels to which convection can penetrate, plus 1.
151     ! NL MUST be less than or equal to ND-1.
152 guez 62
153 guez 52 ! delt: The model time step (sec) between calls to CONVECT
154 guez 62
155 guez 52 !----------------------------------------------------------------------------
156     ! --- On Output:
157 guez 62
158 guez 52 ! iflag: An output integer whose value denotes the following:
159     ! VALUE INTERPRETATION
160     ! ----- --------------
161     ! 0 Moist convection occurs.
162     ! 1 Moist convection occurs, but a CFL condition
163     ! on the subsidence warming is violated. This
164     ! does not cause the scheme to terminate.
165     ! 2 Moist convection, but no precip because ep(inb) lt 0.0001
166     ! 3 No moist convection because new cbmf is 0 and old cbmf is 0.
167     ! 4 No moist convection; atmosphere is not
168     ! unstable
169     ! 6 No moist convection because ihmin le minorig.
170     ! 7 No moist convection because unreasonable
171     ! parcel level temperature or specific humidity.
172     ! 8 No moist convection: lifted condensation
173     ! level is above the 200 mb level.
174     ! 9 No moist convection: cloud base is higher
175     ! then the level NL-1.
176 guez 62
177 guez 52 ! ft: Array of temperature tendency (K/s) of dimension ND, defined at same
178     ! grid levels as T, Q, QS and P.
179 guez 62
180 guez 52 ! fq: Array of specific humidity tendencies ((gm/gm)/s) of dimension ND,
181     ! defined at same grid levels as T, Q, QS and P.
182 guez 62
183 guez 52 ! fu: Array of forcing of zonal velocity (m/s^2) of dimension ND,
184     ! defined at same grid levels as T.
185 guez 62
186 guez 52 ! fv: Same as FU, but for forcing of meridional velocity.
187 guez 62
188 guez 52 ! ftra: Array of forcing of tracer content, in tracer mixing ratio per
189     ! second, defined at same levels as T. Dimensioned (ND, NTRA).
190 guez 62
191 guez 52 ! precip: Scalar convective precipitation rate (mm/day).
192 guez 62
193 guez 52 ! VPrecip: Vertical profile of convective precipitation (kg/m2/s).
194 guez 62
195 guez 52 ! wd: A convective downdraft velocity scale. For use in surface
196     ! flux parameterizations. See convect.ps file for details.
197 guez 62
198 guez 52 ! tprime: A convective downdraft temperature perturbation scale (K).
199     ! For use in surface flux parameterizations. See convect.ps
200     ! file for details.
201 guez 62
202 guez 52 ! qprime: A convective downdraft specific humidity
203     ! perturbation scale (gm/gm).
204     ! For use in surface flux parameterizations. See convect.ps
205     ! file for details.
206 guez 62
207 guez 52 ! cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST
208     ! BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT
209     ! ITS NEXT CALL. That is, the value of CBMF must be "remembered"
210     ! by the calling program between calls to CONVECT.
211 guez 62
212 guez 52 ! det: Array of detrainment mass flux of dimension ND.
213 guez 62
214 guez 52 !-------------------------------------------------------------------
215 guez 62
216 guez 52 ! Local arrays
217 guez 3
218 guez 72 integer noff
219     real da(len, nd), phi(len, nd, nd), mp(len, nd)
220    
221 guez 52 integer i, k, n, il, j
222     integer icbmax
223     integer nk1(klon)
224     integer icbs1(klon)
225 guez 3
226 guez 52 real plcl1(klon)
227     real tnk1(klon)
228     real qnk1(klon)
229     real gznk1(klon)
230     real pnk1(klon)
231     real qsnk1(klon)
232     real pbase1(klon)
233     real buoybase1(klon)
234 guez 3
235 guez 52 real lv1(klon, klev)
236     real cpn1(klon, klev)
237     real tv1(klon, klev)
238     real gz1(klon, klev)
239     real hm1(klon, klev)
240     real h1(klon, klev)
241     real tp1(klon, klev)
242     real tvp1(klon, klev)
243     real clw1(klon, klev)
244     real th1(klon, klev)
245 guez 62
246 guez 52 integer ncum
247 guez 62
248 guez 52 ! (local) compressed fields:
249 guez 62
250 guez 52 integer nloc
251 guez 91 parameter (nloc = klon) ! pour l'instant
252 guez 3
253 guez 52 integer idcum(nloc)
254     integer iflag(nloc), nk(nloc), icb(nloc)
255     integer nent(nloc, klev)
256     integer icbs(nloc)
257     integer inb(nloc), inbis(nloc)
258 guez 3
259 guez 52 real cbmf(nloc), plcl(nloc), tnk(nloc), qnk(nloc), gznk(nloc)
260     real t(nloc, klev), q(nloc, klev), qs(nloc, klev)
261     real u(nloc, klev), v(nloc, klev)
262     real gz(nloc, klev), h(nloc, klev), lv(nloc, klev), cpn(nloc, klev)
263     real p(nloc, klev), ph(nloc, klev+1), tv(nloc, klev), tp(nloc, klev)
264     real clw(nloc, klev)
265     real dph(nloc, klev)
266     real pbase(nloc), buoybase(nloc), th(nloc, klev)
267     real tvp(nloc, klev)
268     real sig(nloc, klev), w0(nloc, klev)
269     real hp(nloc, klev), ep(nloc, klev), sigp(nloc, klev)
270     real frac(nloc), buoy(nloc, klev)
271     real cape(nloc)
272     real m(nloc, klev), ment(nloc, klev, klev), qent(nloc, klev, klev)
273     real uent(nloc, klev, klev), vent(nloc, klev, klev)
274     real ments(nloc, klev, klev), qents(nloc, klev, klev)
275     real sij(nloc, klev, klev), elij(nloc, klev, klev)
276     real qp(nloc, klev), up(nloc, klev), vp(nloc, klev)
277     real wt(nloc, klev), water(nloc, klev), evap(nloc, klev)
278     real b(nloc, klev), ft(nloc, klev), fq(nloc, klev)
279     real fu(nloc, klev), fv(nloc, klev)
280     real upwd(nloc, klev), dnwd(nloc, klev), dnwd0(nloc, klev)
281     real Ma(nloc, klev), mike(nloc, klev), tls(nloc, klev)
282     real tps(nloc, klev), qprime(nloc), tprime(nloc)
283     real precip(nloc)
284     real VPrecip(nloc, klev+1)
285     real tra(nloc, klev, ntra), trap(nloc, klev, ntra)
286     real ftra(nloc, klev, ntra), traent(nloc, klev, klev, ntra)
287     real qcondc(nloc, klev) ! cld
288     real wd(nloc) ! gust
289 guez 3
290 guez 52 !-------------------------------------------------------------------
291     ! --- SET CONSTANTS AND PARAMETERS
292     !-------------------------------------------------------------------
293 guez 3
294 guez 52 ! -- set simulation flags:
295     ! (common cvflag)
296 guez 3
297 guez 52 CALL cv_flag
298    
299     ! -- set thermodynamical constants:
300     ! (common cvthermo)
301    
302 guez 69 CALL cv_thermo
303 guez 52
304     ! -- set convect parameters
305 guez 62
306 guez 52 ! includes microphysical parameters and parameters that
307     ! control the rate of approach to quasi-equilibrium)
308     ! (common cvparam)
309    
310     if (iflag_con.eq.3) then
311     CALL cv3_param(nd, delt)
312     endif
313    
314     if (iflag_con.eq.4) then
315 guez 3 CALL cv_param(nd)
316 guez 52 endif
317 guez 3
318 guez 52 !---------------------------------------------------------------------
319     ! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS
320     !---------------------------------------------------------------------
321 guez 3
322 guez 91 do k = 1, nd
323     do i = 1, len
324     ft1(i, k) = 0.0
325     fq1(i, k) = 0.0
326     fu1(i, k) = 0.0
327     fv1(i, k) = 0.0
328     tvp1(i, k) = 0.0
329     tp1(i, k) = 0.0
330     clw1(i, k) = 0.0
331 guez 52 !ym
332 guez 91 clw(i, k) = 0.0
333     gz1(i, k) = 0.
334 guez 52 VPrecip1(i, k) = 0.
335 guez 91 Ma1(i, k) = 0.0
336     upwd1(i, k) = 0.0
337     dnwd1(i, k) = 0.0
338     dnwd01(i, k) = 0.0
339     qcondc1(i, k) = 0.0
340 guez 52 end do
341     end do
342 guez 3
343 guez 91 do j = 1, ntra
344     do k = 1, nd
345     do i = 1, len
346     ftra1(i, k, j) = 0.0
347 guez 52 end do
348     end do
349     end do
350 guez 3
351 guez 91 do i = 1, len
352     precip1(i) = 0.0
353     iflag1(i) = 0
354     wd1(i) = 0.0
355     cape1(i) = 0.0
356     VPrecip1(i, nd+1) = 0.0
357 guez 52 end do
358 guez 3
359 guez 52 if (iflag_con.eq.3) then
360 guez 91 do il = 1, len
361     sig1(il, nd) = sig1(il, nd) + 1.
362     sig1(il, nd) = min(sig1(il, nd), 12.1)
363 guez 52 enddo
364     endif
365 guez 3
366 guez 52 !--------------------------------------------------------------------
367     ! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
368     !--------------------------------------------------------------------
369 guez 3
370 guez 52 if (iflag_con.eq.3) then
371     CALL cv3_prelim(len, nd, ndp1, t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, &
372     h1, hm1, th1)! nd->na
373     endif
374 guez 3
375 guez 52 if (iflag_con.eq.4) then
376     CALL cv_prelim(len, nd, ndp1, t1, q1, p1, ph1 &
377     , lv1, cpn1, tv1, gz1, h1, hm1)
378     endif
379 guez 3
380 guez 52 !--------------------------------------------------------------------
381     ! --- CONVECTIVE FEED
382     !--------------------------------------------------------------------
383 guez 3
384 guez 52 if (iflag_con.eq.3) then
385     CALL cv3_feed(len, nd, t1, q1, qs1, p1, ph1, hm1, gz1 &
386     , nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! nd->na
387     endif
388 guez 3
389 guez 52 if (iflag_con.eq.4) then
390     CALL cv_feed(len, nd, t1, q1, qs1, p1, hm1, gz1 &
391     , nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1)
392     endif
393 guez 3
394 guez 52 !--------------------------------------------------------------------
395     ! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part
396     ! (up through ICB for convect4, up through ICB+1 for convect3)
397     ! Calculates the lifted parcel virtual temperature at nk, the
398     ! actual temperature, and the adiabatic liquid water content.
399     !--------------------------------------------------------------------
400 guez 3
401 guez 52 if (iflag_con.eq.3) then
402     CALL cv3_undilute1(len, nd, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1 &
403     , tp1, tvp1, clw1, icbs1) ! nd->na
404     endif
405 guez 3
406 guez 52 if (iflag_con.eq.4) then
407     CALL cv_undilute1(len, nd, t1, q1, qs1, gz1, p1, nk1, icb1, icbmax &
408     , tp1, tvp1, clw1)
409     endif
410 guez 3
411 guez 52 !-------------------------------------------------------------------
412     ! --- TRIGGERING
413     !-------------------------------------------------------------------
414 guez 3
415 guez 52 if (iflag_con.eq.3) then
416 guez 72 CALL cv3_trigger(len, nd, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, &
417     buoybase1, iflag1, sig1, w01) ! nd->na
418 guez 52 endif
419 guez 3
420 guez 52 if (iflag_con.eq.4) then
421     CALL cv_trigger(len, nd, icb1, cbmf1, tv1, tvp1, iflag1)
422     endif
423 guez 3
424 guez 52 ! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY
425 guez 3
426 guez 91 ncum = 0
427     do i = 1, len
428 guez 52 if(iflag1(i).eq.0)then
429 guez 91 ncum = ncum+1
430     idcum(ncum) = i
431 guez 52 endif
432     end do
433 guez 3
434 guez 52 ! print*, 'klon, ncum = ', len, ncum
435 guez 3
436 guez 52 IF (ncum.gt.0) THEN
437 guez 3
438 guez 52 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
439     ! --- COMPRESS THE FIELDS
440     ! (-> vectorization over convective gridpoints)
441     !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
442 guez 3
443 guez 52 if (iflag_con.eq.3) then
444 guez 72 CALL cv3_compress(len, nloc, ncum, nd, ntra, iflag1, nk1, icb1, &
445     icbs1, plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1, t1, q1, &
446     qs1, u1, v1, gz1, th1, tra1, h1, lv1, cpn1, p1, ph1, tv1, tp1, &
447     tvp1, clw1, sig1, w01, iflag, nk, icb, icbs, plcl, tnk, qnk, &
448     gznk, pbase, buoybase, t, q, qs, u, v, gz, th, tra, h, lv, &
449     cpn, p, ph, tv, tp, tvp, clw, sig, w0)
450 guez 52 endif
451 guez 3
452 guez 52 if (iflag_con.eq.4) then
453     CALL cv_compress( len, nloc, ncum, nd &
454     , iflag1, nk1, icb1 &
455     , cbmf1, plcl1, tnk1, qnk1, gznk1 &
456     , t1, q1, qs1, u1, v1, gz1 &
457     , h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1 &
458     , iflag, nk, icb &
459     , cbmf, plcl, tnk, qnk, gznk &
460     , t, q, qs, u, v, gz, h, lv, cpn, p, ph, tv, tp, tvp, clw &
461     , dph )
462     endif
463 guez 3
464 guez 52 !-------------------------------------------------------------------
465     ! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :
466     ! --- FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
467     ! --- &
468     ! --- COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
469     ! --- FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
470     ! --- &
471     ! --- FIND THE LEVEL OF NEUTRAL BUOYANCY
472     !-------------------------------------------------------------------
473 guez 3
474 guez 52 if (iflag_con.eq.3) then
475     CALL cv3_undilute2(nloc, ncum, nd, icb, icbs, nk &
476     , tnk, qnk, gznk, t, q, qs, gz &
477     , p, h, tv, lv, pbase, buoybase, plcl &
478     , inb, tp, tvp, clw, hp, ep, sigp, buoy) !na->nd
479     endif
480 guez 3
481 guez 52 if (iflag_con.eq.4) then
482     CALL cv_undilute2(nloc, ncum, nd, icb, nk &
483     , tnk, qnk, gznk, t, q, qs, gz &
484     , p, dph, h, tv, lv &
485     , inb, inbis, tp, tvp, clw, hp, ep, sigp, frac)
486     endif
487 guez 3
488 guez 52 !-------------------------------------------------------------------
489     ! --- CLOSURE
490     !-------------------------------------------------------------------
491 guez 3
492 guez 52 if (iflag_con.eq.3) then
493     CALL cv3_closure(nloc, ncum, nd, icb, inb &
494     , pbase, p, ph, tv, buoy &
495     , sig, w0, cape, m) ! na->nd
496     endif
497 guez 3
498 guez 52 if (iflag_con.eq.4) then
499     CALL cv_closure(nloc, ncum, nd, nk, icb &
500     , tv, tvp, p, ph, dph, plcl, cpn &
501     , iflag, cbmf)
502     endif
503 guez 3
504 guez 52 !-------------------------------------------------------------------
505     ! --- MIXING
506     !-------------------------------------------------------------------
507 guez 3
508 guez 52 if (iflag_con.eq.3) then
509     CALL cv3_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb &
510     , ph, t, q, qs, u, v, tra, h, lv, qnk &
511     , hp, tv, tvp, ep, clw, m, sig &
512     , ment, qent, uent, vent, nent, sij, elij, ments, qents, traent)! na->nd
513     endif
514 guez 3
515 guez 52 if (iflag_con.eq.4) then
516     CALL cv_mixing(nloc, ncum, nd, icb, nk, inb, inbis &
517     , ph, t, q, qs, u, v, h, lv, qnk &
518     , hp, tv, tvp, ep, clw, cbmf &
519     , m, ment, qent, uent, vent, nent, sij, elij)
520     endif
521 guez 3
522 guez 52 !-------------------------------------------------------------------
523     ! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS
524     !-------------------------------------------------------------------
525 guez 3
526 guez 52 if (iflag_con.eq.3) then
527     CALL cv3_unsat(nloc, ncum, nd, nd, ntra, icb, inb &
528     , t, q, qs, gz, u, v, tra, p, ph &
529     , th, tv, lv, cpn, ep, sigp, clw &
530     , m, ment, elij, delt, plcl &
531     , mp, qp, up, vp, trap, wt, water, evap, b)! na->nd
532     endif
533 guez 3
534 guez 52 if (iflag_con.eq.4) then
535     CALL cv_unsat(nloc, ncum, nd, inb, t, q, qs, gz, u, v, p, ph &
536     , h, lv, ep, sigp, clw, m, ment, elij &
537     , iflag, mp, qp, up, vp, wt, water, evap)
538     endif
539 guez 3
540 guez 52 !-------------------------------------------------------------------
541     ! --- YIELD
542     ! (tendencies, precipitation, variables of interface with other
543     ! processes, etc)
544     !-------------------------------------------------------------------
545 guez 3
546 guez 52 if (iflag_con.eq.3) then
547     CALL cv3_yield(nloc, ncum, nd, nd, ntra &
548     , icb, inb, delt &
549     , t, q, u, v, tra, gz, p, ph, h, hp, lv, cpn, th &
550     , ep, clw, m, tp, mp, qp, up, vp, trap &
551     , wt, water, evap, b &
552     , ment, qent, uent, vent, nent, elij, traent, sig &
553     , tv, tvp &
554     , iflag, precip, VPrecip, ft, fq, fu, fv, ftra &
555     , upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc, wd)! na->nd
556     endif
557 guez 3
558 guez 52 if (iflag_con.eq.4) then
559     CALL cv_yield(nloc, ncum, nd, nk, icb, inb, delt &
560     , t, q, u, v, gz, p, ph, h, hp, lv, cpn &
561     , ep, clw, frac, m, mp, qp, up, vp &
562     , wt, water, evap &
563     , ment, qent, uent, vent, nent, elij &
564     , tv, tvp &
565     , iflag, wd, qprime, tprime &
566     , precip, cbmf, ft, fq, fu, fv, Ma, qcondc)
567     endif
568 guez 3
569 guez 52 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
570     ! --- passive tracers
571     !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
572 guez 3
573 guez 52 if (iflag_con.eq.3) then
574     CALL cv3_tracer(nloc, len, ncum, nd, nd, &
575     ment, sij, da, phi)
576     endif
577 guez 3
578 guez 52 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
579     ! --- UNCOMPRESS THE FIELDS
580     !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
581 guez 91 ! set iflag1 = 42 for non convective points
582     do i = 1, len
583     iflag1(i) = 42
584 guez 52 end do
585 guez 62
586 guez 52 if (iflag_con.eq.3) then
587     CALL cv3_uncompress(nloc, len, ncum, nd, ntra, idcum &
588     , iflag &
589     , precip, VPrecip, sig, w0 &
590     , ft, fq, fu, fv, ftra &
591     , inb &
592     , Ma, upwd, dnwd, dnwd0, qcondc, wd, cape &
593     , da, phi, mp &
594     , iflag1 &
595     , precip1, VPrecip1, sig1, w01 &
596     , ft1, fq1, fu1, fv1, ftra1 &
597     , inb1 &
598     , Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, cape1 &
599     , da1, phi1, mp1)
600     endif
601 guez 3
602 guez 52 if (iflag_con.eq.4) then
603     CALL cv_uncompress(nloc, len, ncum, nd, idcum &
604     , iflag &
605     , precip, cbmf &
606     , ft, fq, fu, fv &
607     , Ma, qcondc &
608     , iflag1 &
609     , precip1, cbmf1 &
610     , ft1, fq1, fu1, fv1 &
611     , Ma1, qcondc1 )
612     endif
613     ENDIF ! ncum>0
614 guez 3
615 guez 52 end SUBROUTINE cv_driver
616 guez 3
617 guez 52 end module cv_driver_m

  ViewVC Help
Powered by ViewVC 1.1.21