source: LMDZ4/trunk/libf/phylmd/cva_driver.F @ 1040

Last change on this file since 1040 was 991, checked in by Laurent Fairhead, 16 years ago

La variable nent n'etait pas remontee dans la bonne routine YM/JYG
LF

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