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

  ViewVC Help
Powered by ViewVC 1.1.21