source: LMDZ4/trunk/libf/phylmd/cv_driver.F @ 524

Last change on this file since 524 was 524, checked in by lmdzadmin, 20 years ago

Initial revision

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