source: LMDZ4/branches/LMDZ4_par_0/libf/phylmd/cv_driver.F @ 634

Last change on this file since 634 was 634, checked in by Laurent Fairhead, 19 years ago

Modifications faites à la physique pour la rendre parallele YM
Une branche de travail LMDZ4_par_0 a été créée provisoirement afin de tester
les modifs pleinement avant leurs inclusions dans le tronc principal
LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 23.7 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,VPrecip1,
8     &                   cbmf1,sig1,w01,
9     &                   icb1,inb1,
10     &                   delt,Ma1,upwd1,dnwd1,dnwd01,qcondc1,wd1,cape1,
11     &                   da1,phi1,mp1)
12C
13      USE dimphy
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
42C      VPrecip1      Real           Output       vertical profile of precipitations
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
67cym#include "dimensions.h"
68cym#include "dimphy.h"
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)
90      real VPrecip1(len,nd+1)
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
100      real da1(len,nd),phi1(len,nd,nd),mp1(len,nd)
101      real da(len,nd),phi(len,nd,nd),mp(len,nd)
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!
194!  VPrecip: Vertical profile of convective precipitation (kg/m2/s).
195!
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)
224      integer inb1(klon)
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
253cym      integer nloc
254cym      parameter (nloc=klon) ! pour l'instant
255#define nloc klon
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)
279      real qp(nloc,klev), up(nloc,klev), vp(nloc,klev)
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)
287      real VPrecip(nloc,klev+1)
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
293!-------------------------------------------------------------------
294! --- SET CONSTANTS AND PARAMETERS
295!-------------------------------------------------------------------
296
297c -- set simulation flags:
298c   (common cvflag)
299
300       CALL cv_flag
301
302c -- set thermodynamical constants:
303c       (common cvthermo)
304
305       CALL cv_thermo(iflag_con)
306
307c -- set convect parameters
308c
309c       includes microphysical parameters and parameters that
310c       control the rate of approach to quasi-equilibrium)
311c       (common cvparam)
312
313      if (iflag_con.eq.3) then
314       CALL cv3_param(nd,delt)
315      endif
316
317      if (iflag_con.eq.4) then
318       CALL cv_param(nd)
319      endif
320
321!---------------------------------------------------------------------
322! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS
323!---------------------------------------------------------------------
324
325      do 20 k=1,nd
326        do 10 i=1,len
327         ft1(i,k)=0.0
328         fq1(i,k)=0.0
329         fu1(i,k)=0.0
330         fv1(i,k)=0.0
331         tvp1(i,k)=0.0
332         tp1(i,k)=0.0
333         clw1(i,k)=0.0
334cym
335         clw(i,k)=0.0   
336         gz1(i,k) = 0.
337         VPrecip1(i,k) = 0.
338         Ma1(i,k)=0.0
339         upwd1(i,k)=0.0
340         dnwd1(i,k)=0.0
341         dnwd01(i,k)=0.0
342         qcondc1(i,k)=0.0
343 10     continue
344 20   continue
345
346      do 30 j=1,ntra
347       do 31 k=1,nd
348        do 32 i=1,len
349         ftra1(i,k,j)=0.0
350 32     continue   
351 31    continue   
352 30   continue   
353
354      do 60 i=1,len
355        precip1(i)=0.0
356        iflag1(i)=0
357        wd1(i)=0.0
358        cape1(i)=0.0
359        VPrecip1(i,nd+1)=0.0
360 60   continue
361
362      if (iflag_con.eq.3) then
363        do il=1,len
364         sig1(il,nd)=sig1(il,nd)+1.
365         sig1(il,nd)=amin1(sig1(il,nd),12.1)
366        enddo
367      endif
368
369!--------------------------------------------------------------------
370! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
371!--------------------------------------------------------------------
372
373      if (iflag_con.eq.3) then
374       CALL cv3_prelim(len,nd,ndp1,t1,q1,p1,ph1            ! nd->na
375     o               ,lv1,cpn1,tv1,gz1,h1,hm1,th1)
376      endif
377
378      if (iflag_con.eq.4) then
379       CALL cv_prelim(len,nd,ndp1,t1,q1,p1,ph1
380     o               ,lv1,cpn1,tv1,gz1,h1,hm1)
381      endif
382
383!--------------------------------------------------------------------
384! --- CONVECTIVE FEED
385!--------------------------------------------------------------------
386
387      if (iflag_con.eq.3) then
388       CALL cv3_feed(len,nd,t1,q1,qs1,p1,ph1,hm1,gz1           ! nd->na
389     o         ,nk1,icb1,icbmax,iflag1,tnk1,qnk1,gznk1,plcl1)
390      endif
391
392      if (iflag_con.eq.4) then
393       CALL cv_feed(len,nd,t1,q1,qs1,p1,hm1,gz1
394     o         ,nk1,icb1,icbmax,iflag1,tnk1,qnk1,gznk1,plcl1)
395      endif
396
397!--------------------------------------------------------------------
398! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part
399! (up through ICB for convect4, up through ICB+1 for convect3)
400!     Calculates the lifted parcel virtual temperature at nk, the
401!     actual temperature, and the adiabatic liquid water content.
402!--------------------------------------------------------------------
403
404      if (iflag_con.eq.3) then
405       CALL cv3_undilute1(len,nd,t1,q1,qs1,gz1,plcl1,p1,nk1,icb1  ! nd->na
406     o                        ,tp1,tvp1,clw1,icbs1)
407      endif
408
409      if (iflag_con.eq.4) then
410       CALL cv_undilute1(len,nd,t1,q1,qs1,gz1,p1,nk1,icb1,icbmax
411     :                        ,tp1,tvp1,clw1)
412      endif
413
414!-------------------------------------------------------------------
415! --- TRIGGERING
416!-------------------------------------------------------------------
417
418      if (iflag_con.eq.3) then
419       CALL cv3_trigger(len,nd,icb1,plcl1,p1,th1,tv1,tvp1      ! nd->na
420     o                 ,pbase1,buoybase1,iflag1,sig1,w01)
421      endif
422
423      if (iflag_con.eq.4) then
424       CALL cv_trigger(len,nd,icb1,cbmf1,tv1,tvp1,iflag1)
425      endif
426
427!=====================================================================
428! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY
429!=====================================================================
430
431      ncum=0
432      do 400 i=1,len
433        if(iflag1(i).eq.0)then
434           ncum=ncum+1
435           idcum(ncum)=i
436        endif
437 400  continue
438
439c       print*,'klon, ncum = ',len,ncum
440
441      IF (ncum.gt.0) THEN
442
443!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
444! --- COMPRESS THE FIELDS
445!               (-> vectorization over convective gridpoints)
446!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
447
448      if (iflag_con.eq.3) then
449       CALL cv3_compress( len,nloc,ncum,nd,ntra
450     :    ,iflag1,nk1,icb1,icbs1
451     :    ,plcl1,tnk1,qnk1,gznk1,pbase1,buoybase1
452     :    ,t1,q1,qs1,u1,v1,gz1,th1
453     :    ,tra1
454     :    ,h1,lv1,cpn1,p1,ph1,tv1,tp1,tvp1,clw1
455     :    ,sig1,w01
456     o    ,iflag,nk,icb,icbs
457     o    ,plcl,tnk,qnk,gznk,pbase,buoybase
458     o    ,t,q,qs,u,v,gz,th
459     o    ,tra
460     o    ,h,lv,cpn,p,ph,tv,tp,tvp,clw
461     o    ,sig,w0  )
462      endif
463
464      if (iflag_con.eq.4) then
465       CALL cv_compress( len,nloc,ncum,nd
466     :    ,iflag1,nk1,icb1
467     :    ,cbmf1,plcl1,tnk1,qnk1,gznk1
468     :    ,t1,q1,qs1,u1,v1,gz1
469     :    ,h1,lv1,cpn1,p1,ph1,tv1,tp1,tvp1,clw1
470     o    ,iflag,nk,icb
471     o    ,cbmf,plcl,tnk,qnk,gznk
472     o    ,t,q,qs,u,v,gz,h,lv,cpn,p,ph,tv,tp,tvp,clw
473     o    ,dph )
474      endif
475
476!-------------------------------------------------------------------
477! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :
478! ---   FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
479! ---   &
480! ---   COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
481! ---   FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
482! ---   &
483! ---   FIND THE LEVEL OF NEUTRAL BUOYANCY
484!-------------------------------------------------------------------
485
486      if (iflag_con.eq.3) then
487       CALL cv3_undilute2(nloc,ncum,nd,icb,icbs,nk        !na->nd
488     :                        ,tnk,qnk,gznk,t,q,qs,gz
489     :                        ,p,h,tv,lv,pbase,buoybase,plcl
490     o                        ,inb,tp,tvp,clw,hp,ep,sigp,buoy)
491      endif
492
493      if (iflag_con.eq.4) then
494       CALL cv_undilute2(nloc,ncum,nd,icb,nk
495     :                        ,tnk,qnk,gznk,t,q,qs,gz
496     :                        ,p,dph,h,tv,lv
497     o             ,inb,inbis,tp,tvp,clw,hp,ep,sigp,frac)
498      endif
499
500!-------------------------------------------------------------------
501! --- CLOSURE
502!-------------------------------------------------------------------
503
504      if (iflag_con.eq.3) then
505       CALL cv3_closure(nloc,ncum,nd,icb,inb              ! na->nd
506     :                       ,pbase,p,ph,tv,buoy
507     o                       ,sig,w0,cape,m)
508      endif
509
510      if (iflag_con.eq.4) then
511       CALL cv_closure(nloc,ncum,nd,nk,icb
512     :                ,tv,tvp,p,ph,dph,plcl,cpn
513     o                ,iflag,cbmf)
514      endif
515
516!-------------------------------------------------------------------
517! --- MIXING
518!-------------------------------------------------------------------
519
520      if (iflag_con.eq.3) then
521       CALL cv3_mixing(nloc,ncum,nd,nd,ntra,icb,nk,inb    ! na->nd
522     :                     ,ph,t,q,qs,u,v,tra,h,lv,qnk
523     :                     ,hp,tv,tvp,ep,clw,m,sig
524     o ,ment,qent,uent,vent,sij,elij,ments,qents,traent)
525      endif
526
527      if (iflag_con.eq.4) then
528       CALL cv_mixing(nloc,ncum,nd,icb,nk,inb,inbis
529     :                     ,ph,t,q,qs,u,v,h,lv,qnk
530     :                     ,hp,tv,tvp,ep,clw,cbmf
531     o                     ,m,ment,qent,uent,vent,nent,sij,elij)
532      endif
533
534!-------------------------------------------------------------------
535! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS
536!-------------------------------------------------------------------
537
538      if (iflag_con.eq.3) then
539       CALL cv3_unsat(nloc,ncum,nd,nd,ntra,icb,inb    ! na->nd
540     :               ,t,q,qs,gz,u,v,tra,p,ph
541     :               ,th,tv,lv,cpn,ep,sigp,clw
542     :               ,m,ment,elij,delt,plcl
543     o          ,mp,qp,up,vp,trap,wt,water,evap,b)
544      endif
545
546      if (iflag_con.eq.4) then
547       CALL cv_unsat(nloc,ncum,nd,inb,t,q,qs,gz,u,v,p,ph
548     :                   ,h,lv,ep,sigp,clw,m,ment,elij
549     o                   ,iflag,mp,qp,up,vp,wt,water,evap)
550      endif
551
552!-------------------------------------------------------------------
553! --- YIELD
554!     (tendencies, precipitation, variables of interface with other
555!      processes, etc)
556!-------------------------------------------------------------------
557
558      if (iflag_con.eq.3) then
559       CALL cv3_yield(nloc,ncum,nd,nd,ntra            ! na->nd
560     :                     ,icb,inb,delt
561     :                     ,t,q,u,v,tra,gz,p,ph,h,hp,lv,cpn,th
562     :                     ,ep,clw,m,tp,mp,qp,up,vp,trap
563     :                     ,wt,water,evap,b
564     :                     ,ment,qent,uent,vent,nent,elij,traent,sig
565     :                     ,tv,tvp
566     o                     ,iflag,precip,VPrecip,ft,fq,fu,fv,ftra
567     o                     ,upwd,dnwd,dnwd0,ma,mike,tls,tps,qcondc,wd)
568      endif
569
570      if (iflag_con.eq.4) then
571       CALL cv_yield(nloc,ncum,nd,nk,icb,inb,delt
572     :              ,t,q,u,v,gz,p,ph,h,hp,lv,cpn
573     :              ,ep,clw,frac,m,mp,qp,up,vp
574     :              ,wt,water,evap
575     :              ,ment,qent,uent,vent,nent,elij
576     :              ,tv,tvp
577     o              ,iflag,wd,qprime,tprime
578     o              ,precip,cbmf,ft,fq,fu,fv,Ma,qcondc)
579      endif
580
581!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
582! --- passive tracers
583!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
584
585      if (iflag_con.eq.3) then
586       CALL cv3_tracer(nloc,len,ncum,nd,nd,
587     :                  ment,sij,da,phi)
588      endif
589
590!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
591! --- UNCOMPRESS THE FIELDS
592!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
593c set iflag1 =42 for non convective points
594      do  i=1,len
595        iflag1(i)=42
596      end do
597c
598      if (iflag_con.eq.3) then
599       CALL cv3_uncompress(nloc,len,ncum,nd,ntra,idcum
600     :          ,iflag
601     :          ,precip,VPrecip,sig,w0
602     :          ,ft,fq,fu,fv,ftra
603     :          ,inb
604     :          ,Ma,upwd,dnwd,dnwd0,qcondc,wd,cape
605     :          ,da,phi,mp
606     o          ,iflag1
607     o          ,precip1,VPrecip1,sig1,w01
608     o          ,ft1,fq1,fu1,fv1,ftra1
609     o          ,inb1
610     o          ,Ma1,upwd1,dnwd1,dnwd01,qcondc1,wd1,cape1
611     o          ,da1,phi1,mp1)
612      endif
613
614      if (iflag_con.eq.4) then
615       CALL cv_uncompress(nloc,len,ncum,nd,idcum
616     :          ,iflag
617     :          ,precip,cbmf
618     :          ,ft,fq,fu,fv
619     :          ,Ma,qcondc           
620     o          ,iflag1
621     o          ,precip1,cbmf1
622     o          ,ft1,fq1,fu1,fv1
623     o          ,Ma1,qcondc1 )           
624      endif
625
626      ENDIF ! ncum>0
627
6289999  continue
629
630      return
631      end
632
633!==================================================================
634      SUBROUTINE cv_flag
635      implicit none
636
637#include "cvflag.h"
638
639c -- si .TRUE., on rend la gravite plus explicite et eventuellement
640c differente de 10.0 dans convect3:
641      cvflag_grav = .TRUE.
642
643      return
644      end
645
646!==================================================================
647      SUBROUTINE cv_thermo(iflag_con)
648          implicit none
649
650c-------------------------------------------------------------
651c Set thermodynamical constants for convectL
652c-------------------------------------------------------------
653
654#include "YOMCST.h"
655#include "cvthermo.h"
656
657      integer iflag_con
658
659
660c original set from convect:
661      if (iflag_con.eq.4) then
662       cpd=1005.7
663       cpv=1870.0
664       cl=4190.0
665       rrv=461.5
666       rrd=287.04
667       lv0=2.501E6
668       g=9.8
669       t0=273.15
670       grav=g
671      endif
672
673c constants consistent with LMDZ:
674      if (iflag_con.eq.3) then
675       cpd = RCPD
676       cpv = RCPV
677       cl  = RCW
678       rrv = RV
679       rrd = RD
680       lv0 = RLVTT
681       g   = RG     ! not used in convect3
682c ori      t0  = RTT
683       t0  = 273.15 ! convect3 (RTT=273.16)
684c maf       grav= 10.    ! implicitely or explicitely used in convect3
685       grav= g    ! implicitely or explicitely used in convect3
686      endif
687
688      rowl=1000.0 !(a quelle variable de YOMCST cela correspond-il?)
689
690      clmcpv=cl-cpv
691      clmcpd=cl-cpd
692      cpdmcp=cpd-cpv
693      cpvmcpd=cpv-cpd
694      cpvmcl=cl-cpv ! for convect3
695      eps=rrd/rrv
696      epsi=1.0/eps
697      epsim1=epsi-1.0
698c      ginv=1.0/g
699      ginv=1.0/grav
700      hrd=0.5*rrd
701
702      return
703      end
704
Note: See TracBrowser for help on using the repository browser.