/[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 69 - (hide annotations)
Mon Feb 18 16:33:12 2013 UTC (11 years, 3 months ago) by guez
Original Path: trunk/libf/phylmd/cv_driver.f90
File size: 23615 byte(s)
Deleted files cvparam3.f90 and nuagecom.f90. Moved variables from
module cvparam3 to module cv3_param_m. Moved variables rad_chau1 and
rad_chau2 from module nuagecom to module conf_phys_m.

Read clesphys2_nml from conf_phys instead of gcm.

Removed argument iflag_con from several procedures. Access module
variable instead.

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

  ViewVC Help
Powered by ViewVC 1.1.21