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

Last change on this file since 5360 was 965, checked in by Laurent Fairhead, 17 years ago

Un peu de menage pour avoir moins de sorties a l'execution
LF

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