source: LMDZ4/branches/LMDZ4-dev-20091210/libf/phylmd/cva_driver.F @ 1440

Last change on this file since 1440 was 1248, checked in by Laurent Fairhead, 15 years ago

Correction de bug JYG

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 34.1 KB
Line 
1      SUBROUTINE cva_driver(len,nd,ndp1,ntra,nloc,
2     &                   iflag_con,iflag_mix,
3     &                   iflag_clos,delt,
4     &                   t1,q1,qs1,t1_wake,q1_wake,qs1_wake,s1_wake,
5     &                   u1,v1,tra1,
6     &                   p1,ph1,
7     &                   ALE1,ALP1,
8     &                   sig1feed1,sig2feed1,wght1,
9     o                   iflag1,ft1,fq1,fu1,fv1,ftra1,
10     &                   precip1,kbas1,ktop1,cbmf1,
11     &                   sig1,w01,                  !input/output
12     &                   ptop21,sigd,
13     &                   Ma1,mip1,Vprecip1,upwd1,dnwd1,dnwd01,
14     &                   qcondc1,wd1,
15     &                   cape1,cin1,tvp1,
16     &                   ftd1,fqd1,
17     &                   Plim11,Plim21,asupmax1,supmax01,asupmaxmin1
18     &                   ,lalim_conv)
19***************************************************************
20*                                                             *
21* CV_DRIVER                                                   *
22*                                                             *
23*                                                             *
24* written by   : Sandrine Bony-Lena , 17/05/2003, 11.19.41    *
25* modified by :                                               *
26***************************************************************
27***************************************************************
28C
29      USE dimphy
30      implicit none
31C
32C.............................START PROLOGUE............................
33C
34C PARAMETERS:
35C      Name            Type         Usage            Description
36C   ----------      ----------     -------  ----------------------------
37C
38C      len           Integer        Input        first (i) dimension
39C      nd            Integer        Input        vertical (k) dimension
40C      ndp1          Integer        Input        nd + 1
41C      ntra          Integer        Input        number of tracors
42C      iflag_con     Integer        Input        version of convect (3/4)
43C      iflag_mix     Integer        Input        version of mixing  (0/1/2)
44C      iflag_clos    Integer        Input        version of closure (0/1)
45C      delt          Real           Input        time step
46C      t1            Real           Input        temperature (sat draught envt)
47C      q1            Real           Input        specific hum (sat draught envt)
48C      qs1           Real           Input        sat specific hum (sat draught envt)
49C      t1_wake       Real           Input        temperature (unsat draught envt)
50C      q1_wake       Real           Input        specific hum(unsat draught envt)
51C      qs1_wake      Real           Input        sat specific hum(unsat draughts envt)
52C      s1_wake       Real           Input        fractionnal area covered by wakes
53C      u1            Real           Input        u-wind
54C      v1            Real           Input        v-wind
55C      tra1          Real           Input        tracors
56C      p1            Real           Input        full level pressure
57C      ph1           Real           Input        half level pressure
58C      ALE1          Real           Input        Available lifting Energy
59C      ALP1          Real           Input        Available lifting Power
60C      sig1feed1     Real           Input        sigma coord at lower bound of feeding layer
61C      sig2feed1     Real           Input        sigma coord at upper bound of feeding layer
62C      wght1         Real           Input        weight density determining the feeding mixture
63C      iflag1        Integer        Output       flag for Emanuel conditions
64C      ft1           Real           Output       temp tend
65C      fq1           Real           Output       spec hum tend
66C      fu1           Real           Output       u-wind tend
67C      fv1           Real           Output       v-wind tend
68C      ftra1         Real           Output       tracor tend
69C      precip1       Real           Output       precipitation
70C      kbas1         Integer        Output       cloud base level
71C      ktop1         Integer        Output       cloud top level
72C      cbmf1         Real           Output       cloud base mass flux
73C      sig1          Real           In/Out       section adiabatic updraft
74C      w01           Real           In/Out       vertical velocity within adiab updraft
75C      ptop21        Real           In/Out       top of entraining zone
76C      Ma1           Real           Output       mass flux adiabatic updraft
77C      mip1          Real           Output       mass flux shed by the adiabatic updraft
78C      Vprecip1      Real           Output       vertical profile of precipitations
79C      upwd1         Real           Output       total upward mass flux (adiab+mixed)
80C      dnwd1         Real           Output       saturated downward mass flux (mixed)
81C      dnwd01        Real           Output       unsaturated downward mass flux
82C      qcondc1       Real           Output       in-cld mixing ratio of condensed water
83C      wd1           Real           Output       downdraft velocity scale for sfc fluxes
84C      cape1         Real           Output       CAPE
85C      cin1          Real           Output       CIN
86C      tvp1          Real           Output       adiab lifted parcell virt temp
87C      ftd1          Real           Output       precip temp tend
88C      fqt1          Real           Output       precip spec hum tend
89C      Plim11        Real           Output
90C      Plim21        Real           Output
91C      asupmax1      Real           Output
92C      supmax01      Real           Output
93C      asupmaxmin1   Real           Output
94C S. Bony, Mar 2002:
95C       * Several modules corresponding to different physical processes
96C       * Several versions of convect may be used:
97C               - iflag_con=3: version lmd  (previously named convect3)
98C               - iflag_con=4: version 4.3b (vect. version, previously convect1/2)
99C   + tard:     - iflag_con=5: version lmd with ice (previously named convectg)
100C S. Bony, Oct 2002:
101C       * Vectorization of convect3 (ie version lmd)
102C
103C..............................END PROLOGUE.............................
104c
105c
106#include "dimensions.h"
107ccccc#include "dimphy.h"
108c
109c Input
110      integer len
111      integer nd
112      integer ndp1
113      integer ntra
114      integer iflag_con
115      integer iflag_mix
116      integer iflag_clos
117      real delt
118      real t1(len,nd)
119      real q1(len,nd)
120      real qs1(len,nd)
121      real t1_wake(len,nd)
122      real q1_wake(len,nd)
123      real qs1_wake(len,nd)
124      real s1_wake(len)
125      real u1(len,nd)
126      real v1(len,nd)
127      real tra1(len,nd,ntra)
128      real p1(len,nd)
129      real ph1(len,ndp1)
130      real ALE1(len)
131      real ALP1(len)
132      real sig1feed1 ! pressure at lower bound of feeding layer
133      real sig2feed1 ! pressure at upper bound of feeding layer
134      real wght1(nd) ! weight density determining the feeding mixture
135c
136c Output
137      integer iflag1(len)
138      real ft1(len,nd)
139      real fq1(len,nd)
140      real fu1(len,nd)
141      real fv1(len,nd)
142      real ftra1(len,nd,ntra)
143      real precip1(len)
144      integer kbas1(len)
145      integer ktop1(len)
146      real cbmf1(len)
147!      real cbmflast(len)
148      real sig1(len,klev)      !input/output
149      real w01(len,klev)       !input/output
150      real ptop21(len)
151      real Ma1(len,nd)
152      real mip1(len,nd)
153      real Vprecip1(len,nd)
154      real upwd1(len,nd)
155      real dnwd1(len,nd)
156      real dnwd01(len,nd)
157      real qcondc1(len,nd)     ! cld
158      real wd1(len)            ! gust
159      real cape1(len)
160      real cin1(len)
161      real tvp1(len,nd)
162c
163      real ftd1(len,nd)
164      real fqd1(len,nd)
165      real Plim11(len)
166      real Plim21(len)
167      real asupmax1(len,nd)
168      real supmax01(len)
169      real asupmaxmin1(len)
170      integer lalim_conv(len)
171!-------------------------------------------------------------------
172! --- ARGUMENTS
173!-------------------------------------------------------------------
174! --- On input:
175!
176!  t:   Array of absolute temperature (K) of dimension ND, with first
177!       index corresponding to lowest model level. Note that this array
178!       will be altered by the subroutine if dry convective adjustment
179!       occurs and if IPBL is not equal to 0.
180!
181!  q:   Array of specific humidity (gm/gm) of dimension ND, with first
182!       index corresponding to lowest model level. Must be defined
183!       at same grid levels as T. Note that this array will be altered
184!       if dry convective adjustment occurs and if IPBL is not equal to 0.
185!
186!  qs:  Array of saturation specific humidity of dimension ND, with first
187!       index corresponding to lowest model level. Must be defined
188!       at same grid levels as T. Note that this array will be altered
189!       if dry convective adjustment occurs and if IPBL is not equal to 0.
190!
191! t_wake: Array of absolute temperature (K), seen by unsaturated draughts,
192!       of dimension ND, with first index corresponding to lowest model level.
193!
194! q_wake: Array of specific humidity (gm/gm), seen by unsaturated draughts,
195!       of dimension ND, with first index corresponding to lowest model level.
196!       Must be defined at same grid levels as T.
197!
198!qs_wake: Array of saturation specific humidity, seen by unsaturated draughts,
199!       of dimension ND, with first index corresponding to lowest model level.
200!       Must be defined at same grid levels as T.
201!
202!s_wake: Array of fractionnal area occupied by the wakes.
203!
204!  u:   Array of zonal wind velocity (m/s) of dimension ND, witth first
205!       index corresponding with the lowest model level. Defined at
206!       same levels as T. Note that this array will be altered if
207!       dry convective adjustment occurs and if IPBL is not equal to 0.
208!
209!  v:   Same as u but for meridional velocity.
210!
211!  tra: Array of passive tracer mixing ratio, of dimensions (ND,NTRA),
212!       where NTRA is the number of different tracers. If no
213!       convective tracer transport is needed, define a dummy
214!       input array of dimension (ND,1). Tracers are defined at
215!       same vertical levels as T. Note that this array will be altered
216!       if dry convective adjustment occurs and if IPBL is not equal to 0.
217!
218!  p:   Array of pressure (mb) of dimension ND, with first
219!       index corresponding to lowest model level. Must be defined
220!       at same grid levels as T.
221!
222!  ph:  Array of pressure (mb) of dimension ND+1, with first index
223!       corresponding to lowest level. These pressures are defined at
224!       levels intermediate between those of P, T, Q and QS. The first
225!       value of PH should be greater than (i.e. at a lower level than)
226!       the first value of the array P.
227!
228! ALE:  Available lifting Energy
229!
230! ALP:  Available lifting Power
231!
232!  nl:  The maximum number of levels to which convection can penetrate, plus 1.
233!       NL MUST be less than or equal to ND-1.
234!
235!  delt: The model time step (sec) between calls to CONVECT
236!
237!----------------------------------------------------------------------------
238! ---   On Output:
239!
240!  iflag: An output integer whose value denotes the following:
241!       VALUE   INTERPRETATION
242!       -----   --------------
243!         0     Moist convection occurs.
244!         1     Moist convection occurs, but a CFL condition
245!               on the subsidence warming is violated. This
246!               does not cause the scheme to terminate.
247!         2     Moist convection, but no precip because ep(inb) lt 0.0001
248!         3     No moist convection because new cbmf is 0 and old cbmf is 0.
249!         4     No moist convection; atmosphere is not
250!               unstable
251!         6     No moist convection because ihmin le minorig.
252!         7     No moist convection because unreasonable
253!               parcel level temperature or specific humidity.
254!         8     No moist convection: lifted condensation
255!               level is above the 200 mb level.
256!         9     No moist convection: cloud base is higher
257!               then the level NL-1.
258!
259!  ft:   Array of temperature tendency (K/s) of dimension ND, defined at same
260!        grid levels as T, Q, QS and P.
261!
262!  fq:   Array of specific humidity tendencies ((gm/gm)/s) of dimension ND,
263!        defined at same grid levels as T, Q, QS and P.
264!
265!  fu:   Array of forcing of zonal velocity (m/s^2) of dimension ND,
266!        defined at same grid levels as T.
267!
268!  fv:   Same as FU, but for forcing of meridional velocity.
269!
270!  ftra: Array of forcing of tracer content, in tracer mixing ratio per
271!        second, defined at same levels as T. Dimensioned (ND,NTRA).
272!
273!  precip: Scalar convective precipitation rate (mm/day).
274!
275!  wd:   A convective downdraft velocity scale. For use in surface
276!        flux parameterizations. See convect.ps file for details.
277!
278!  tprime: A convective downdraft temperature perturbation scale (K).
279!          For use in surface flux parameterizations. See convect.ps
280!          file for details.
281!
282!  qprime: A convective downdraft specific humidity
283!          perturbation scale (gm/gm).
284!          For use in surface flux parameterizations. See convect.ps
285!          file for details.
286!
287!  cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST
288!        BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT
289!        ITS NEXT CALL. That is, the value of CBMF must be "remembered"
290!        by the calling program between calls to CONVECT.
291!
292!  det:   Array of detrainment mass flux of dimension ND.
293!
294!  ftd:  Array of temperature tendency due to precipitations (K/s) of dimension ND,
295!        defined at same grid levels as T, Q, QS and P.
296!
297!  fqd:  Array of specific humidity tendencies due to precipitations ((gm/gm)/s)
298!        of dimension ND, defined at same grid levels as T, Q, QS and P.
299!
300!-------------------------------------------------------------------
301c
302c  Local arrays
303c
304
305      integer i,k,n,il,j
306      integer nword1,nword2,nword3,nword4
307      integer icbmax
308      integer nk1(klon)
309      integer icb1(klon)
310      integer icbs1(klon)
311
312      logical ok_inhib  ! True => possible inhibition of convection by dryness
313      logical, save :: debut=.true.
314c$OMP THREADPRIVATE(debut)
315
316      real plcl1(klon)
317      real tnk1(klon)
318      real thnk1(klon)
319      real qnk1(klon)
320      real gznk1(klon)
321      real pnk1(klon)
322      real qsnk1(klon)
323      real unk1(klon)
324      real vnk1(klon)
325      real cpnk1(klon)
326      real hnk1(klon)
327      real pbase1(klon)
328      real buoybase1(klon)
329
330      real lv1(klon,klev) ,lv1_wake(klon,klev)
331      real cpn1(klon,klev),cpn1_wake(klon,klev)
332      real tv1(klon,klev) ,tv1_wake(klon,klev)
333      real gz1(klon,klev) ,gz1_wake(klon,klev)
334      real hm1(klon,klev) ,hm1_wake(klon,klev)
335      real h1(klon,klev)  ,h1_wake(klon,klev)
336      real tp1(klon,klev)
337      real clw1(klon,klev)
338      real th1(klon,klev) ,th1_wake(klon,klev)
339c
340      real bid(klon,klev)   ! dummy array
341c
342      integer ncum
343c
344      integer j1feed(klon)
345      integer j2feed(klon)
346      real p1feed1(len) ! pressure at lower bound of feeding layer
347      real p2feed1(len) ! pressure at upper bound of feeding layer
348      real wghti1(len,nd) ! weights of the feeding layers
349c
350c (local) compressed fields:
351c
352      integer nloc
353c      parameter (nloc=klon) ! pour l'instant
354
355      integer idcum(nloc)
356      integer iflag(nloc),nk(nloc),icb(nloc)
357      integer nent(nloc,klev)
358      integer icbs(nloc)
359      integer inb(nloc), inbis(nloc)
360
361      real cbmf(nloc),plcl(nloc)
362      real t(nloc,klev),q(nloc,klev),qs(nloc,klev)
363      real t_wake(nloc,klev),q_wake(nloc,klev),qs_wake(nloc,klev)
364      real s_wake(nloc)
365      real u(nloc,klev),v(nloc,klev)
366      real gz(nloc,klev),h(nloc,klev)     ,hm(nloc,klev)
367      real               h_wake(nloc,klev),hm_wake(nloc,klev)
368      real lv(nloc,klev)     ,cpn(nloc,klev)
369      real lv_wake(nloc,klev),cpn_wake(nloc,klev)
370      real p(nloc,klev),ph(nloc,klev+1),tv(nloc,klev)    ,tp(nloc,klev)
371      real                              tv_wake(nloc,klev)
372      real clw(nloc,klev)
373      real dph(nloc,klev)
374      real pbase(nloc), buoybase(nloc), th(nloc,klev)
375      real                              th_wake(nloc,klev)
376      real tvp(nloc,klev)
377      real sig(nloc,klev), w0(nloc,klev), ptop2(nloc)
378      real hp(nloc,klev), ep(nloc,klev), sigp(nloc,klev)
379      real frac(nloc), buoy(nloc,klev)
380      real cape(nloc)
381      real cin(nloc)
382      real m(nloc,klev)
383      real ment(nloc,klev,klev), sij(nloc,klev,klev)
384      real qent(nloc,klev,klev)
385      real hent(nloc,klev,klev)
386      real uent(nloc,klev,klev), vent(nloc,klev,klev)
387      real ments(nloc,klev,klev), qents(nloc,klev,klev)
388      real elij(nloc,klev,klev)
389      real supmax(nloc,klev)
390      real ale(nloc),alp(nloc),coef_clos(nloc)
391      real sigd(nloc)
392!      real mp(nloc,klev), qp(nloc,klev), up(nloc,klev), vp(nloc,klev)
393!      real wt(nloc,klev), water(nloc,klev), evap(nloc,klev)
394!      real b(nloc,klev), sigd(nloc)
395!      save mp,qp,up,vp,wt,water,evap,b
396      real, save, allocatable :: mp(:,:),qp(:,:),up(:,:),vp(:,:)
397      real, save, allocatable :: wt(:,:),water(:,:),evap(:,:), b(:,:)
398c$OMP THREADPRIVATE(mp,qp,up,vp,wt,water,evap,b)
399      real  ft(nloc,klev), fq(nloc,klev)
400      real ftd(nloc,klev), fqd(nloc,klev)
401      real fu(nloc,klev), fv(nloc,klev)
402      real upwd(nloc,klev), dnwd(nloc,klev), dnwd0(nloc,klev)
403      real Ma(nloc,klev), mip(nloc,klev), tls(nloc,klev)
404      real tps(nloc,klev), qprime(nloc), tprime(nloc)
405      real precip(nloc)
406      real Vprecip(nloc,klev)
407      real tra(nloc,klev,ntra), trap(nloc,klev,ntra)
408      real ftra(nloc,klev,ntra), traent(nloc,klev,klev,ntra)
409      real qcondc(nloc,klev)  ! cld
410      real wd(nloc)           ! gust
411      real Plim1(nloc),Plim2(nloc)
412      real asupmax(nloc,klev)
413      real supmax0(nloc)
414      real asupmaxmin(nloc)
415c
416      real tnk(nloc),qnk(nloc),gznk(nloc)
417      real wghti(nloc,nd)
418      real hnk(nloc),unk(nloc),vnk(nloc)
419      logical, save :: first=.true.
420c$OMP THREADPRIVATE(first)
421
422c
423!      print *, 't1, t1_wake ',(k,t1(1,k),t1_wake(1,k),k=1,klev)
424!      print *, 'q1, q1_wake ',(k,q1(1,k),q1_wake(1,k),k=1,klev)
425
426!-------------------------------------------------------------------
427! --- SET CONSTANTS AND PARAMETERS
428!-------------------------------------------------------------------
429
430       if (first) then
431         allocate(mp(nloc,klev), qp(nloc,klev), up(nloc,klev))
432         allocate(vp(nloc,klev), wt(nloc,klev), water(nloc,klev))
433         allocate(evap(nloc,klev), b(nloc,klev))
434         first=.false.
435       endif
436c -- set simulation flags:
437c   (common cvflag)
438
439       CALL cv_flag
440
441c -- set thermodynamical constants:
442c       (common cvthermo)
443
444       CALL cv_thermo(iflag_con)
445
446c -- set convect parameters
447c
448c       includes microphysical parameters and parameters that
449c       control the rate of approach to quasi-equilibrium)
450c       (common cvparam)
451
452      if (iflag_con.eq.3) then
453       CALL cv3_param(nd,delt)
454 
455      endif
456
457      if (iflag_con.eq.4) then
458       CALL cv_param(nd)
459      endif
460
461!---------------------------------------------------------------------
462! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS
463!---------------------------------------------------------------------
464      nword1=len
465      nword2=len*nd
466      nword3=len*nd*ntra
467      nword4=len*nd*nd
468 
469!      call izilch(iflag1  ,nword1)
470!      call  zilch(iflag1  ,nword1)
471      do i=1,len
472         iflag1(i)=0
473         ktop1(i)=0
474         kbas1(i)=0
475      enddo
476      call  zilch(ft1     ,nword2)
477      call  zilch(fq1     ,nword2)
478      call  zilch(fu1     ,nword2)
479      call  zilch(fv1     ,nword2)
480      call  zilch(ftra1   ,nword3)
481      call  zilch(precip1 ,nword1)
482!      call izilch(kbas1   ,nword1)
483!      call  zilch(kbas1   ,nword1)
484!      call izilch(ktop1   ,nword1)
485!      call  zilch(ktop1   ,nword1)
486      call  zilch(cbmf1   ,nword1)
487      call  zilch(ptop21  ,nword1)
488      call  zilch(Ma1     ,nword2)
489      call  zilch(mip1    ,nword2)
490      call  zilch(Vprecip1,nword2)
491      call  zilch(upwd1   ,nword2)
492      call  zilch(dnwd1   ,nword2)
493      call  zilch(dnwd01  ,nword2)
494      call  zilch(qcondc1 ,nword2)
495!test
496!      call  zilch(qcondc ,nword2)
497      call  zilch(wd1     ,nword1)
498      call  zilch(cape1   ,nword1)
499      call  zilch(cin1    ,nword1)
500      call  zilch(tvp1    ,nword2)
501      call  zilch(ftd1    ,nword2)
502      call  zilch(fqd1    ,nword2)
503      call  zilch(Plim11  ,nword1)
504      call  zilch(Plim21  ,nword1)
505      call  zilch(asupmax1,nword2)
506      call  zilch(supmax01,nword1)
507      call  zilch(asupmaxmin1,nword1)
508c
509      DO il = 1,len
510       cin1(il) = -100000.
511       cape1(il) = -1.
512      ENDDO
513
514      if (iflag_con.eq.3) then
515        do il=1,len
516         sig1(il,nd)=sig1(il,nd)+1.
517         sig1(il,nd)=amin1(sig1(il,nd),12.1)
518        enddo
519      endif
520 
521!---------------------------------------------------------------------
522! --- INITIALIZE LOCAL ARRAYS AND PARAMETERS
523!---------------------------------------------------------------------
524!
525      do il = 1,nloc
526         coef_clos(il)=1.
527      enddo
528
529!--------------------------------------------------------------------
530! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
531!--------------------------------------------------------------------
532
533      if (iflag_con.eq.3) then
534 
535       if (debut) THEN
536        print*,'Emanuel version 3 nouvelle'
537       endif
538!       print*,'t1, q1 ',t1,q1
539       CALL cv3_prelim(len,nd,ndp1,t1,q1,p1,ph1      ! nd->na
540     o               ,lv1,cpn1,tv1,gz1,h1,hm1,th1)
541   
542c
543       CALL cv3_prelim(len,nd,ndp1,t1_wake,q1_wake,p1,ph1 ! nd->na
544     o               ,lv1_wake,cpn1_wake,tv1_wake,gz1_wake
545     o               ,h1_wake,bid,th1_wake)
546   
547      endif
548c
549      if (iflag_con.eq.4) then
550       print*,'Emanuel version 4 '
551       CALL cv_prelim(len,nd,ndp1,t1,q1,p1,ph1
552     o               ,lv1,cpn1,tv1,gz1,h1,hm1)
553      endif
554
555!--------------------------------------------------------------------
556! --- CONVECTIVE FEED
557!--------------------------------------------------------------------
558!
559! compute feeding layer potential temperature and mixing ratio :
560!
561! get bounds of feeding layer
562!
563c test niveaux couche alimentation KE
564       if(sig1feed1.eq.sig2feed1) then
565               print*,'impossible de choisir sig1feed=sig2feed'
566               print*,'changer la valeur de sig2feed dans physiq.def'
567       stop
568       endif
569c
570       do i=1,len
571         p1feed1(i)=sig1feed1*ph1(i,1)
572         p2feed1(i)=sig2feed1*ph1(i,1)
573ctest maf
574c         p1feed1(i)=ph1(i,1)
575c         p2feed1(i)=ph1(i,2)
576c         p2feed1(i)=ph1(i,3)
577ctestCR: on prend la couche alim des thermiques
578c          p2feed1(i)=ph1(i,lalim_conv(i)+1)
579c          print*,'lentr=',lentr(i),ph1(i,lentr(i)+1),ph1(i,2)
580       end do
581!
582       if (iflag_con.eq.3) then
583       endif
584      do i=1,len
585!      print*,'avant cv3_feed plim',p1feed1(i),p2feed1(i)
586      enddo
587      if (iflag_con.eq.3) then
588 
589c     print*, 'IFLAG1 avant cv3_feed'
590c     print*,'len,nd',len,nd
591c     write(*,'(64i1)') iflag1(2:klon-1)
592
593       CALL cv3_feed(len,nd,t1,q1,u1,v1,p1,ph1,hm1,gz1           ! nd->na
594     i         ,p1feed1,p2feed1,wght1
595     o         ,wghti1,tnk1,thnk1,qnk1,qsnk1,unk1,vnk1
596     o         ,cpnk1,hnk1,nk1,icb1,icbmax,iflag1,gznk1,plcl1)
597      endif
598   
599c     print*, 'IFLAG1 apres cv3_feed'
600c     print*,'len,nd',len,nd
601c     write(*,'(64i1)') iflag1(2:klon-1)
602
603      if (iflag_con.eq.4) then
604       CALL cv_feed(len,nd,t1,q1,qs1,p1,hm1,gz1
605     o         ,nk1,icb1,icbmax,iflag1,tnk1,qnk1,gznk1,plcl1)
606      endif
607c
608!      print *, 'cv3_feed-> iflag1, plcl1 ',iflag1(1),plcl1(1)
609c
610!--------------------------------------------------------------------
611! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part
612! (up through ICB for convect4, up through ICB+1 for convect3)
613!     Calculates the lifted parcel virtual temperature at nk, the
614!     actual temperature, and the adiabatic liquid water content.
615!--------------------------------------------------------------------
616
617      if (iflag_con.eq.3) then
618   
619       CALL cv3_undilute1(len,nd,t1,qs1,gz1,plcl1,p1,icb1,tnk1,qnk1  ! nd->na
620     o                    ,gznk1,tp1,tvp1,clw1,icbs1)
621      endif
622   
623
624      if (iflag_con.eq.4) then
625       CALL cv_undilute1(len,nd,t1,q1,qs1,gz1,p1,nk1,icb1,icbmax
626     :                        ,tp1,tvp1,clw1)
627      endif
628c
629!-------------------------------------------------------------------
630! --- TRIGGERING
631!-------------------------------------------------------------------
632c
633!      print *,' avant triggering, iflag_con ',iflag_con
634c
635      if (iflag_con.eq.3) then
636   
637       CALL cv3_trigger(len,nd,icb1,plcl1,p1,th1,tv1,tvp1,thnk1      ! nd->na
638     o                 ,pbase1,buoybase1,iflag1,sig1,w01)
639   
640
641c     print*, 'IFLAG1 apres cv3_triger'
642c     print*,'len,nd',len,nd
643c     write(*,'(64i1)') iflag1(2:klon-1)
644
645c     call dump2d(iim,jjm-1,sig1(2)
646      endif
647
648      if (iflag_con.eq.4) then
649       CALL cv_trigger(len,nd,icb1,cbmf1,tv1,tvp1,iflag1)
650      endif
651c
652c
653!=====================================================================
654! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY
655!=====================================================================
656
657      ncum=0
658      do 400 i=1,len
659        if(iflag1(i).eq.0)then
660           ncum=ncum+1
661           idcum(ncum)=i
662        endif
663 400  continue
664c
665!       print*,'klon, ncum = ',len,ncum
666c
667      IF (ncum.gt.0) THEN
668
669!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
670! --- COMPRESS THE FIELDS
671!               (-> vectorization over convective gridpoints)
672!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
673
674      if (iflag_con.eq.3) then
675!       print*,'ncum tv1 ',ncum,tv1
676!       print*,'tvp1 ',tvp1
677       CALL cv3a_compress( len,nloc,ncum,nd,ntra
678     :    ,iflag1,nk1,icb1,icbs1
679     :    ,plcl1,tnk1,qnk1,gznk1,hnk1,unk1,vnk1
680     :    ,wghti1,pbase1,buoybase1
681     :    ,t1,q1,qs1,t1_wake,q1_wake,qs1_wake,s1_wake
682     :    ,u1,v1,gz1,th1,th1_wake
683     :    ,tra1
684     :    ,h1     ,lv1     ,cpn1   ,p1,ph1,tv1    ,tp1,tvp1,clw1
685     :    ,h1_wake,lv1_wake,cpn1_wake     ,tv1_wake
686     :    ,sig1,w01,ptop21
687     :    ,Ale1,Alp1
688     o    ,iflag,nk,icb,icbs
689     o    ,plcl,tnk,qnk,gznk,hnk,unk,vnk
690     o    ,wghti,pbase,buoybase
691     o    ,t,q,qs,t_wake,q_wake,qs_wake,s_wake
692     o    ,u,v,gz,th,th_wake
693     o    ,tra
694     o    ,h     ,lv     ,cpn    ,p,ph,tv    ,tp,tvp,clw
695     o    ,h_wake,lv_wake,cpn_wake    ,tv_wake
696     o    ,sig,w0,ptop2
697     o    ,Ale,Alp  )
698
699!       print*,'tv ',tv
700!       print*,'tvp ',tvp
701
702      endif
703
704      if (iflag_con.eq.4) then
705       CALL cv_compress( len,nloc,ncum,nd
706     :    ,iflag1,nk1,icb1
707     :    ,cbmf1,plcl1,tnk1,qnk1,gznk1
708     :    ,t1,q1,qs1,u1,v1,gz1
709     :    ,h1,lv1,cpn1,p1,ph1,tv1,tp1,tvp1,clw1
710     o    ,iflag,nk,icb
711     o    ,cbmf,plcl,tnk,qnk,gznk
712     o    ,t,q,qs,u,v,gz,h,lv,cpn,p,ph,tv,tp,tvp,clw
713     o    ,dph )
714      endif
715
716!-------------------------------------------------------------------
717! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :
718! ---   FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
719! ---   &
720! ---   COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
721! ---   FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
722! ---   &
723! ---   FIND THE LEVEL OF NEUTRAL BUOYANCY
724!-------------------------------------------------------------------
725
726      if (iflag_con.eq.3) then
727       CALL cv3_undilute2(nloc,ncum,nd,icb,icbs,nk        !na->nd
728     :                        ,tnk,qnk,gznk,hnk,t,q,qs,gz
729     :                        ,p,h,tv,lv,pbase,buoybase,plcl
730     o                        ,inb,tp,tvp,clw,hp,ep,sigp,buoy)
731   
732      endif
733
734      if (iflag_con.eq.4) then
735       CALL cv_undilute2(nloc,ncum,nd,icb,nk
736     :                        ,tnk,qnk,gznk,t,q,qs,gz
737     :                        ,p,dph,h,tv,lv
738     o             ,inb,inbis,tp,tvp,clw,hp,ep,sigp,frac)
739      endif
740c
741!-------------------------------------------------------------------
742! --- MIXING(1)   (if iflag_mix .ge. 1)
743!-------------------------------------------------------------------
744      IF (iflag_con .eq. 3) THEN
745       IF (iflag_mix .ge. 1 ) THEN
746         CALL zilch(supmax,nloc*klev)   
747         CALL cv3p_mixing(nloc,ncum,nd,nd,ntra,icb,nk,inb    ! na->nd
748     :                       ,ph,t,q,qs,u,v,tra,h,lv,qnk
749     :                       ,unk,vnk,hp,tv,tvp,ep,clw,sig
750     :                    ,ment,qent,hent,uent,vent,nent
751     :                   ,sij,elij,supmax,ments,qents,traent)
752!        print*, 'cv3p_mixing-> supmax ', (supmax(1,k), k=1,nd)
753     
754       ELSE
755        CALL zilch(supmax,nloc*klev)
756       ENDIF
757      ENDIF
758!-------------------------------------------------------------------
759! --- CLOSURE
760!-------------------------------------------------------------------
761
762c
763      if (iflag_con.eq.3) then
764       IF (iflag_clos .eq. 0) THEN
765        CALL cv3_closure(nloc,ncum,nd,icb,inb              ! na->nd
766     :                       ,pbase,p,ph,tv,buoy
767     o                       ,sig,w0,cape,m,iflag)
768       ENDIF
769c
770       ok_inhib = iflag_mix .EQ. 2
771c
772       IF (iflag_clos .eq. 1) THEN
773        print *,' pas d appel cv3p_closure'
774cc        CALL cv3p_closure(nloc,ncum,nd,icb,inb              ! na->nd
775cc    :                       ,pbase,plcl,p,ph,tv,tvp,buoy
776cc    :                       ,supmax
777cc    o                       ,sig,w0,ptop2,cape,cin,m)
778       ENDIF
779       IF (iflag_clos .eq. 2) THEN
780        CALL cv3p1_closure(nloc,ncum,nd,icb,inb              ! na->nd
781     :                       ,pbase,plcl,p,ph,tv,tvp,buoy
782     :                       ,supmax,ok_inhib,Ale,Alp
783     o                       ,sig,w0,ptop2,cape,cin,m,iflag,coef_clos
784     :                       ,Plim1,Plim2,asupmax,supmax0
785     :                       ,asupmaxmin,cbmf)
786       ENDIF
787      endif   ! iflag_con.eq.3
788 
789      if (iflag_con.eq.4) then
790       CALL cv_closure(nloc,ncum,nd,nk,icb
791     :                ,tv,tvp,p,ph,dph,plcl,cpn
792     o                ,iflag,cbmf)
793      endif
794c
795!      print *,'cv_closure-> cape ',cape(1)
796c
797!-------------------------------------------------------------------
798! --- MIXING(2)
799!-------------------------------------------------------------------
800
801      if (iflag_con.eq.3) then
802        IF (iflag_mix.eq.0) THEN
803         CALL cv3_mixing(nloc,ncum,nd,nd,ntra,icb,nk,inb    ! na->nd
804     :                       ,ph,t,q,qs,u,v,tra,h,lv,qnk
805     :                       ,unk,vnk,hp,tv,tvp,ep,clw,m,sig
806     o   ,ment,qent,uent,vent,nent,sij,elij,ments,qents,traent)
807         CALL zilch(hent,nloc*klev*klev)
808        ELSE
809         CALL cv3_mixscale(nloc,ncum,nd,ment,m)
810         if (debut) THEN
811          print *,' cv3_mixscale-> '
812         endif !(debut) THEN
813        ENDIF
814      endif
815
816      if (iflag_con.eq.4) then
817       CALL cv_mixing(nloc,ncum,nd,icb,nk,inb,inbis
818     :                     ,ph,t,q,qs,u,v,h,lv,qnk
819     :                     ,hp,tv,tvp,ep,clw,cbmf
820     o                     ,m,ment,qent,uent,vent,nent,sij,elij)
821      endif
822c
823      if (debut) THEN
824       print *,' cv_mixing ->'
825      endif !(debut) THEN
826c      do i = 1,klev
827c        print*,'cv_mixing-> i,ment ',i,(ment(1,i,j),j=1,klev)
828c      enddo
829c
830!-------------------------------------------------------------------
831! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS
832!-------------------------------------------------------------------
833      if (iflag_con.eq.3) then
834       if (debut) THEN
835        print *,' cva_driver -> cv3_unsat '
836       endif !(debut) THEN
837   
838       CALL cv3_unsat(nloc,ncum,nd,nd,ntra,icb,inb,iflag    ! na->nd
839     :               ,t_wake,q_wake,qs_wake,gz,u,v,tra,p,ph
840     :               ,th_wake,tv_wake,lv_wake,cpn_wake
841     :               ,ep,sigp,clw
842     :               ,m,ment,elij,delt,plcl,coef_clos
843     o          ,mp,qp,up,vp,trap,wt,water,evap,b,sigd)
844      endif
845     
846      if (iflag_con.eq.4) then
847       CALL cv_unsat(nloc,ncum,nd,inb,t,q,qs,gz,u,v,p,ph
848     :                   ,h,lv,ep,sigp,clw,m,ment,elij
849     o                   ,iflag,mp,qp,up,vp,wt,water,evap)
850      endif
851c
852      if (debut) THEN
853       print *,'cv_unsat-> '
854       debut=.FALSE.
855      endif !(debut) THEN
856!
857c      print *,'cv_unsat-> mp ',mp
858c      print *,'cv_unsat-> water ',water
859!-------------------------------------------------------------------
860! --- YIELD
861!     (tendencies, precipitation, variables of interface with other
862!      processes, etc)
863!-------------------------------------------------------------------
864
865      if (iflag_con.eq.3) then
866 
867       CALL cv3_yield(nloc,ncum,nd,nd,ntra            ! na->nd
868     :                     ,icb,inb,delt
869     :                     ,t,q,t_wake,q_wake,s_wake,u,v,tra
870     :                     ,gz,p,ph,h,hp,lv,cpn,th,th_wake
871     :                     ,ep,clw,m,tp,mp,qp,up,vp,trap
872     :                     ,wt,water,evap,b,sigd
873     :                    ,ment,qent,hent,iflag_mix,uent,vent
874     :                    ,nent,elij,traent,sig
875     :                    ,tv,tvp,wghti
876     :                    ,iflag,precip,Vprecip,ft,fq,fu,fv,ftra
877     :                    ,cbmf,upwd,dnwd,dnwd0,ma,mip
878     :                    ,tls,tps,qcondc,wd
879     :                    ,ftd,fqd)
880!      print *,' cv3_yield -> fqd(1) = ',fqd(1,1)
881      endif
882
883      if (iflag_con.eq.4) then
884       CALL cv_yield(nloc,ncum,nd,nk,icb,inb,delt
885     :              ,t,q,t_wake,q_wake,u,v,tra
886     :              ,gz,p,ph,h,hp,lv,cpn,th
887     :              ,ep,clw,frac,m,mp,qp,up,vp
888     :              ,wt,water,evap
889     :              ,ment,qent,uent,vent,nent,elij
890     :              ,tv,tvp
891     o              ,iflag,wd,qprime,tprime
892     o              ,precip,cbmf,ft,fq,fu,fv,Ma,qcondc)
893      endif
894
895!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
896! --- UNCOMPRESS THE FIELDS
897!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
898
899
900      if (iflag_con.eq.3) then
901       CALL cv3a_uncompress(nloc,len,ncum,nd,ntra,idcum
902     :          ,iflag,icb,inb
903     :          ,precip,sig,w0,ptop2
904     :          ,ft,fq,fu,fv,ftra
905     :          ,Ma,mip,Vprecip,upwd,dnwd,dnwd0
906     ;          ,qcondc,wd,cape,cin
907     :          ,tvp
908     :          ,ftd,fqd
909     :          ,Plim1,Plim2,asupmax,supmax0
910     :          ,asupmaxmin
911     o          ,iflag1,kbas1,ktop1
912     o          ,precip1,sig1,w01,ptop21
913     o          ,ft1,fq1,fu1,fv1,ftra1
914     o          ,Ma1,mip1,Vprecip1,upwd1,dnwd1,dnwd01
915     o          ,qcondc1,wd1,cape1,cin1
916     o          ,tvp1
917     o          ,ftd1,fqd1
918     o          ,Plim11,Plim21,asupmax1,supmax01
919     o          ,asupmaxmin1)
920      endif
921
922      if (iflag_con.eq.4) then
923       CALL cv_uncompress(nloc,len,ncum,nd,idcum
924     :          ,iflag
925     :          ,precip,cbmf
926     :          ,ft,fq,fu,fv
927     :          ,Ma,qcondc
928     o          ,iflag1
929     o          ,precip1,cbmf1
930     o          ,ft1,fq1,fu1,fv1
931     o          ,Ma1,qcondc1 )
932      endif
933
934      ENDIF ! ncum>0
935
9369999  continue
937      return
938      end
939
Note: See TracBrowser for help on using the repository browser.