source: LMDZ.3.3/branches/rel-LF/libf/phylmd/cv_driver.F @ 1546

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

Suite des adaptations au portage sur SGI et VPP pour Prism (Caubel & Demory)
LF

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