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

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

On remplace le fichier include dimphy.h par le module dimphy.F90i pour etre
coherent avec le partout
LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 33.6 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,
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      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"
106ccccc#include "dimphy.h"
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
309      logical, save :: debut=.true.
310
311      real plcl1(klon)
312      real tnk1(klon)
313      real thnk1(klon)
314      real qnk1(klon)
315      real gznk1(klon)
316      real pnk1(klon)
317      real qsnk1(klon)
318      real unk1(klon)
319      real vnk1(klon)
320      real cpnk1(klon)
321      real hnk1(klon)
322      real pbase1(klon)
323      real buoybase1(klon)
324
325      real lv1(klon,klev) ,lv1_wake(klon,klev)
326      real cpn1(klon,klev),cpn1_wake(klon,klev)
327      real tv1(klon,klev) ,tv1_wake(klon,klev)
328      real gz1(klon,klev) ,gz1_wake(klon,klev)
329      real hm1(klon,klev) ,hm1_wake(klon,klev)
330      real h1(klon,klev)  ,h1_wake(klon,klev)
331      real tp1(klon,klev)
332      real clw1(klon,klev)
333      real th1(klon,klev) ,th1_wake(klon,klev)
334c
335      real bid(klon,klev)   ! dummy array
336c
337      integer ncum
338c
339      integer j1feed(klon)
340      integer j2feed(klon)
341      real p1feed1(len) ! pressure at lower bound of feeding layer
342      real p2feed1(len) ! pressure at upper bound of feeding layer
343      real wghti1(len,nd) ! weights of the feeding layers
344c
345c (local) compressed fields:
346c
347      integer nloc
348c      parameter (nloc=klon) ! pour l'instant
349
350      integer idcum(nloc)
351      integer iflag(nloc),nk(nloc),icb(nloc)
352      integer nent(nloc,klev)
353      integer icbs(nloc)
354      integer inb(nloc), inbis(nloc)
355
356      real cbmf(nloc),plcl(nloc)
357      real t(nloc,klev),q(nloc,klev),qs(nloc,klev)
358      real t_wake(nloc,klev),q_wake(nloc,klev),qs_wake(nloc,klev)
359      real u(nloc,klev),v(nloc,klev)
360      real gz(nloc,klev),h(nloc,klev)     ,hm(nloc,klev)
361      real               h_wake(nloc,klev),hm_wake(nloc,klev)
362      real lv(nloc,klev)     ,cpn(nloc,klev)
363      real lv_wake(nloc,klev),cpn_wake(nloc,klev)
364      real p(nloc,klev),ph(nloc,klev+1),tv(nloc,klev)    ,tp(nloc,klev)
365      real                              tv_wake(nloc,klev)
366      real clw(nloc,klev)
367      real dph(nloc,klev)
368      real pbase(nloc), buoybase(nloc), th(nloc,klev)
369      real                              th_wake(nloc,klev)
370      real tvp(nloc,klev)
371      real sig(nloc,klev), w0(nloc,klev), ptop2(nloc)
372      real hp(nloc,klev), ep(nloc,klev), sigp(nloc,klev)
373      real frac(nloc), buoy(nloc,klev)
374      real cape(nloc)
375      real cin(nloc)
376      real m(nloc,klev)
377      real ment(nloc,klev,klev), sij(nloc,klev,klev)
378      real qent(nloc,klev,klev)
379      real hent(nloc,klev,klev)
380      real uent(nloc,klev,klev), vent(nloc,klev,klev)
381      real ments(nloc,klev,klev), qents(nloc,klev,klev)
382      real elij(nloc,klev,klev)
383      real supmax(nloc,klev)
384      real ale(nloc),alp(nloc),coef_clos(nloc)
385      real sigd(nloc)
386!      real mp(nloc,klev), qp(nloc,klev), up(nloc,klev), vp(nloc,klev)
387!      real wt(nloc,klev), water(nloc,klev), evap(nloc,klev)
388!      real b(nloc,klev), sigd(nloc)
389!      save mp,qp,up,vp,wt,water,evap,b
390      real, save, allocatable :: mp(:,:),qp(:,:),up(:,:),vp(:,:)
391      real, save, allocatable :: wt(:,:),water(:,:),evap(:,:), b(:,:)
392c$OMP THREADPRIVATE(mp,qp,up,vp,wt,water,evap,b)
393      real  ft(nloc,klev), fq(nloc,klev)
394      real ftd(nloc,klev), fqd(nloc,klev)
395      real fu(nloc,klev), fv(nloc,klev)
396      real upwd(nloc,klev), dnwd(nloc,klev), dnwd0(nloc,klev)
397      real Ma(nloc,klev), mip(nloc,klev), tls(nloc,klev)
398      real tps(nloc,klev), qprime(nloc), tprime(nloc)
399      real precip(nloc)
400      real Vprecip(nloc,klev)
401      real tra(nloc,klev,ntra), trap(nloc,klev,ntra)
402      real ftra(nloc,klev,ntra), traent(nloc,klev,klev,ntra)
403      real qcondc(nloc,klev)  ! cld
404      real wd(nloc)           ! gust
405      real Plim1(nloc),Plim2(nloc)
406      real asupmax(nloc,klev)
407      real supmax0(nloc)
408      real asupmaxmin(nloc)
409c
410      real tnk(nloc),qnk(nloc),gznk(nloc)
411      real wghti(nloc,nd)
412      real hnk(nloc),unk(nloc),vnk(nloc)
413      logical, save :: first=.true.
414
415c
416!      print *, 't1, t1_wake ',(k,t1(1,k),t1_wake(1,k),k=1,klev)
417!      print *, 'q1, q1_wake ',(k,q1(1,k),q1_wake(1,k),k=1,klev)
418
419!-------------------------------------------------------------------
420! --- SET CONSTANTS AND PARAMETERS
421!-------------------------------------------------------------------
422
423       if (first) then
424         allocate(mp(nloc,klev), qp(nloc,klev), up(nloc,klev))
425         allocate(vp(nloc,klev), wt(nloc,klev), water(nloc,klev))
426         allocate(evap(nloc,klev), b(nloc,klev))
427         first=.false.
428       endif
429c -- set simulation flags:
430c   (common cvflag)
431
432       CALL cv_flag
433
434c -- set thermodynamical constants:
435c       (common cvthermo)
436
437       CALL cv_thermo(iflag_con)
438
439c -- set convect parameters
440c
441c       includes microphysical parameters and parameters that
442c       control the rate of approach to quasi-equilibrium)
443c       (common cvparam)
444
445      if (iflag_con.eq.3) then
446       CALL cv3_param(nd,delt)
447 
448      endif
449
450      if (iflag_con.eq.4) then
451       CALL cv_param(nd)
452      endif
453
454!---------------------------------------------------------------------
455! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS
456!---------------------------------------------------------------------
457      nword1=len
458      nword2=len*nd
459      nword3=len*nd*ntra
460      nword4=len*nd*nd
461 
462!      call izilch(iflag1  ,nword1)
463!      call  zilch(iflag1  ,nword1)
464      do i=1,len
465         iflag1(i)=0
466         ktop1(i)=0
467         kbas1(i)=0
468      enddo
469      call  zilch(ft1     ,nword2)
470      call  zilch(fq1     ,nword2)
471      call  zilch(fu1     ,nword2)
472      call  zilch(fv1     ,nword2)
473      call  zilch(ftra1   ,nword3)
474      call  zilch(precip1 ,nword1)
475!      call izilch(kbas1   ,nword1)
476!      call  zilch(kbas1   ,nword1)
477!      call izilch(ktop1   ,nword1)
478!      call  zilch(ktop1   ,nword1)
479      call  zilch(cbmf1   ,nword1)
480      call  zilch(ptop21  ,nword1)
481      call  zilch(Ma1     ,nword2)
482      call  zilch(mip1    ,nword2)
483      call  zilch(Vprecip1,nword2)
484      call  zilch(upwd1   ,nword2)
485      call  zilch(dnwd1   ,nword2)
486      call  zilch(dnwd01  ,nword2)
487      call  zilch(qcondc1 ,nword2)
488!test
489!      call  zilch(qcondc ,nword2)
490      call  zilch(wd1     ,nword1)
491      call  zilch(cape1   ,nword1)
492      call  zilch(cin1    ,nword1)
493      call  zilch(tvp1    ,nword2)
494      call  zilch(ftd1    ,nword2)
495      call  zilch(fqd1    ,nword2)
496      call  zilch(Plim11  ,nword1)
497      call  zilch(Plim21  ,nword1)
498      call  zilch(asupmax1,nword2)
499      call  zilch(supmax01,nword1)
500      call  zilch(asupmaxmin1,nword1)
501c
502 
503      if (iflag_con.eq.3) then
504        do il=1,len
505         sig1(il,nd)=sig1(il,nd)+1.
506         sig1(il,nd)=amin1(sig1(il,nd),12.1)
507        enddo
508      endif
509 
510!---------------------------------------------------------------------
511! --- INITIALIZE LOCAL ARRAYS AND PARAMETERS
512!---------------------------------------------------------------------
513!
514      do il = 1,nloc
515         coef_clos(il)=1.
516      enddo
517
518!--------------------------------------------------------------------
519! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
520!--------------------------------------------------------------------
521
522      if (iflag_con.eq.3) then
523 
524       if (debut) THEN
525        print*,'Emanuel version 3 nouvelle'
526       endif
527
528       CALL cv3_prelim(len,nd,ndp1,t1,q1,p1,ph1      ! nd->na
529     o               ,lv1,cpn1,tv1,gz1,h1,hm1,th1)
530   
531c
532       CALL cv3_prelim(len,nd,ndp1,t1_wake,q1_wake,p1,ph1 ! nd->na
533     o               ,lv1_wake,cpn1_wake,tv1_wake,gz1_wake
534     o               ,h1_wake,bid,th1_wake)
535   
536      endif
537c
538      if (iflag_con.eq.4) then
539       print*,'Emanuel version 4 '
540       CALL cv_prelim(len,nd,ndp1,t1,q1,p1,ph1
541     o               ,lv1,cpn1,tv1,gz1,h1,hm1)
542      endif
543
544!--------------------------------------------------------------------
545! --- CONVECTIVE FEED
546!--------------------------------------------------------------------
547!
548! compute feeding layer potential temperature and mixing ratio :
549!
550! get bounds of feeding layer
551!
552c test niveaux couche alimentation KE
553       if(sig1feed1.eq.sig2feed1) then
554               print*,'impossible de choisir sig1feed=sig2feed'
555               print*,'changer la valeur de sig2feed dans physiq.def'
556       stop
557       endif
558c
559       do i=1,len
560         p1feed1(i)=sig1feed1*ph1(i,1)
561         p2feed1(i)=sig2feed1*ph1(i,1)
562ctest maf
563c         p1feed1(i)=ph1(i,1)
564c         p2feed1(i)=ph1(i,2)
565c         p2feed1(i)=ph1(i,3)
566ctestCR: on prend la couche alim des thermiques
567c          p2feed1(i)=ph1(i,lalim_conv(i)+1)
568c          print*,'lentr=',lentr(i),ph1(i,lentr(i)+1),ph1(i,2)
569       end do
570!
571       if (iflag_con.eq.3) then
572       endif
573      do i=1,len
574!      print*,'avant cv3_feed plim',p1feed1(i),p2feed1(i)
575      enddo
576      if (iflag_con.eq.3) then
577 
578c     print*, 'IFLAG1 avant cv3_feed'
579c     print*,'len,nd',len,nd
580c     write(*,'(64i1)') iflag1(2:klon-1)
581
582       CALL cv3_feed(len,nd,t1,q1,u1,v1,p1,ph1,hm1,gz1           ! nd->na
583     i         ,p1feed1,p2feed1,wght1
584     o         ,wghti1,tnk1,thnk1,qnk1,qsnk1,unk1,vnk1
585     o         ,cpnk1,hnk1,nk1,icb1,icbmax,iflag1,gznk1,plcl1)
586      endif
587   
588c     print*, 'IFLAG1 apres cv3_feed'
589c     print*,'len,nd',len,nd
590c     write(*,'(64i1)') iflag1(2:klon-1)
591
592      if (iflag_con.eq.4) then
593       CALL cv_feed(len,nd,t1,q1,qs1,p1,hm1,gz1
594     o         ,nk1,icb1,icbmax,iflag1,tnk1,qnk1,gznk1,plcl1)
595      endif
596c
597!      print *, 'cv3_feed-> iflag1, plcl1 ',iflag1(1),plcl1(1)
598c
599!--------------------------------------------------------------------
600! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part
601! (up through ICB for convect4, up through ICB+1 for convect3)
602!     Calculates the lifted parcel virtual temperature at nk, the
603!     actual temperature, and the adiabatic liquid water content.
604!--------------------------------------------------------------------
605
606      if (iflag_con.eq.3) then
607   
608       CALL cv3_undilute1(len,nd,t1,qs1,gz1,plcl1,p1,icb1,tnk1,qnk1  ! nd->na
609     o                    ,gznk1,tp1,tvp1,clw1,icbs1)
610      endif
611   
612
613      if (iflag_con.eq.4) then
614       CALL cv_undilute1(len,nd,t1,q1,qs1,gz1,p1,nk1,icb1,icbmax
615     :                        ,tp1,tvp1,clw1)
616      endif
617c
618!-------------------------------------------------------------------
619! --- TRIGGERING
620!-------------------------------------------------------------------
621c
622!      print *,' avant triggering, iflag_con ',iflag_con
623c
624      if (iflag_con.eq.3) then
625   
626       CALL cv3_trigger(len,nd,icb1,plcl1,p1,th1,tv1,tvp1,thnk1      ! nd->na
627     o                 ,pbase1,buoybase1,iflag1,sig1,w01)
628   
629
630c     print*, 'IFLAG1 apres cv3_triger'
631c     print*,'len,nd',len,nd
632c     write(*,'(64i1)') iflag1(2:klon-1)
633
634c     call dump2d(iim,jjm-1,sig1(2)
635      endif
636
637      if (iflag_con.eq.4) then
638       CALL cv_trigger(len,nd,icb1,cbmf1,tv1,tvp1,iflag1)
639      endif
640c
641c
642!=====================================================================
643! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY
644!=====================================================================
645
646      ncum=0
647      do 400 i=1,len
648        if(iflag1(i).eq.0)then
649           ncum=ncum+1
650           idcum(ncum)=i
651        endif
652 400  continue
653c
654!       print*,'klon, ncum = ',len,ncum
655c
656      IF (ncum.gt.0) THEN
657
658!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
659! --- COMPRESS THE FIELDS
660!               (-> vectorization over convective gridpoints)
661!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
662
663      if (iflag_con.eq.3) then
664     
665       CALL cv3a_compress( len,nloc,ncum,nd,ntra
666     :    ,iflag1,nk1,icb1,icbs1
667     :    ,plcl1,tnk1,qnk1,gznk1,hnk1,unk1,vnk1
668     :    ,wghti1,pbase1,buoybase1
669     :    ,t1,q1,qs1,t1_wake,q1_wake,qs1_wake,u1,v1,gz1,th1,th1_wake
670     :    ,tra1
671     :    ,h1     ,lv1     ,cpn1   ,p1,ph1,tv1    ,tp1,tvp1,clw1
672     :    ,h1_wake,lv1_wake,cpn1_wake     ,tv1_wake
673     :    ,sig1,w01,ptop21
674     :    ,Ale1,Alp1
675     o    ,iflag,nk,icb,icbs
676     o    ,plcl,tnk,qnk,gznk,hnk,unk,vnk
677     o    ,wghti,pbase,buoybase
678     o    ,t,q,qs,t_wake,q_wake,qs_wake,u,v,gz,th,th_wake
679     o    ,tra
680     o    ,h     ,lv     ,cpn    ,p,ph,tv    ,tp,tvp,clw
681     o    ,h_wake,lv_wake,cpn_wake    ,tv_wake
682     o    ,sig,w0,ptop2
683     o    ,Ale,Alp  )
684
685      endif
686
687      if (iflag_con.eq.4) then
688       CALL cv_compress( len,nloc,ncum,nd
689     :    ,iflag1,nk1,icb1
690     :    ,cbmf1,plcl1,tnk1,qnk1,gznk1
691     :    ,t1,q1,qs1,u1,v1,gz1
692     :    ,h1,lv1,cpn1,p1,ph1,tv1,tp1,tvp1,clw1
693     o    ,iflag,nk,icb
694     o    ,cbmf,plcl,tnk,qnk,gznk
695     o    ,t,q,qs,u,v,gz,h,lv,cpn,p,ph,tv,tp,tvp,clw
696     o    ,dph )
697      endif
698
699!-------------------------------------------------------------------
700! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :
701! ---   FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
702! ---   &
703! ---   COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
704! ---   FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
705! ---   &
706! ---   FIND THE LEVEL OF NEUTRAL BUOYANCY
707!-------------------------------------------------------------------
708
709      if (iflag_con.eq.3) then
710       CALL cv3_undilute2(nloc,ncum,nd,icb,icbs,nk        !na->nd
711     :                        ,tnk,qnk,gznk,hnk,t,q,qs,gz
712     :                        ,p,h,tv,lv,pbase,buoybase,plcl
713     o                        ,inb,tp,tvp,clw,hp,ep,sigp,buoy)
714   
715      endif
716
717      if (iflag_con.eq.4) then
718       CALL cv_undilute2(nloc,ncum,nd,icb,nk
719     :                        ,tnk,qnk,gznk,t,q,qs,gz
720     :                        ,p,dph,h,tv,lv
721     o             ,inb,inbis,tp,tvp,clw,hp,ep,sigp,frac)
722      endif
723c
724!-------------------------------------------------------------------
725! --- MIXING(1)   (if iflag_mix .ge. 1)
726!-------------------------------------------------------------------
727      IF (iflag_con .eq. 3) THEN
728       IF (iflag_mix .ge. 1 ) THEN
729   
730         CALL cv3p_mixing(nloc,ncum,nd,nd,ntra,icb,nk,inb    ! na->nd
731     :                       ,ph,t,q,qs,u,v,tra,h,lv,qnk
732     :                       ,unk,vnk,hp,tv,tvp,ep,clw,sig
733     :                    ,ment,qent,hent,uent,vent
734     :                   ,sij,elij,supmax,ments,qents,traent)
735!        print*, 'cv3p_mixing-> supmax ', (supmax(1,k), k=1,nd)
736     
737       ELSE
738        CALL zilch(supmax,nloc*klev)
739       ENDIF
740      ENDIF
741!-------------------------------------------------------------------
742! --- CLOSURE
743!-------------------------------------------------------------------
744
745c
746      if (iflag_con.eq.3) then
747       IF (iflag_clos .eq. 0) THEN
748        CALL cv3_closure(nloc,ncum,nd,icb,inb              ! na->nd
749     :                       ,pbase,p,ph,tv,buoy
750     o                       ,sig,w0,cape,m,iflag)
751       ENDIF
752c
753       ok_inhib = iflag_mix .EQ. 2
754c
755       IF (iflag_clos .eq. 1) THEN
756        print *,' pas d appel cv3p_closure'
757cc        CALL cv3p_closure(nloc,ncum,nd,icb,inb              ! na->nd
758cc    :                       ,pbase,plcl,p,ph,tv,tvp,buoy
759cc    :                       ,supmax
760cc    o                       ,sig,w0,ptop2,cape,cin,m)
761       ENDIF
762       IF (iflag_clos .eq. 2) THEN
763        CALL cv3p1_closure(nloc,ncum,nd,icb,inb              ! na->nd
764     :                       ,pbase,plcl,p,ph,tv,tvp,buoy
765     :                       ,supmax,ok_inhib,Ale,Alp
766     o                       ,sig,w0,ptop2,cape,cin,m,iflag,coef_clos
767     :                       ,Plim1,Plim2,asupmax,supmax0
768     :                       ,asupmaxmin,cbmf1)
769       ENDIF
770      endif   ! iflag_con.eq.3
771 
772      if (iflag_con.eq.4) then
773       CALL cv_closure(nloc,ncum,nd,nk,icb
774     :                ,tv,tvp,p,ph,dph,plcl,cpn
775     o                ,iflag,cbmf)
776      endif
777c
778!      print *,'cv_closure-> cape ',cape(1)
779c
780!-------------------------------------------------------------------
781! --- MIXING(2)
782!-------------------------------------------------------------------
783
784      if (iflag_con.eq.3) then
785        IF (iflag_mix.eq.0) THEN
786         CALL cv3_mixing(nloc,ncum,nd,nd,ntra,icb,nk,inb    ! na->nd
787     :                       ,ph,t,q,qs,u,v,tra,h,lv,qnk
788     :                       ,unk,vnk,hp,tv,tvp,ep,clw,m,sig
789     o   ,ment,qent,uent,vent,sij,elij,ments,qents,traent)
790         CALL zilch(hent,nloc*klev*klev)
791        ELSE
792         CALL cv3_mixscale(nloc,ncum,nd,ment,m)
793         if (debut) THEN
794          print *,' cv3_mixscale-> '
795         endif !(debut) THEN
796        ENDIF
797      endif
798
799      if (iflag_con.eq.4) then
800       CALL cv_mixing(nloc,ncum,nd,icb,nk,inb,inbis
801     :                     ,ph,t,q,qs,u,v,h,lv,qnk
802     :                     ,hp,tv,tvp,ep,clw,cbmf
803     o                     ,m,ment,qent,uent,vent,nent,sij,elij)
804      endif
805c
806      if (debut) THEN
807       print *,' cv_mixing ->'
808      endif !(debut) THEN
809c      do i = 1,klev
810c        print*,'cv_mixing-> i,ment ',i,(ment(1,i,j),j=1,klev)
811c      enddo
812c
813!-------------------------------------------------------------------
814! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS
815!-------------------------------------------------------------------
816      if (iflag_con.eq.3) then
817       if (debut) THEN
818        print *,' cva_driver -> cv3_unsat '
819       endif !(debut) THEN
820   
821       CALL cv3_unsat(nloc,ncum,nd,nd,ntra,icb,inb,iflag    ! na->nd
822     :               ,t_wake,q_wake,qs_wake,gz,u,v,tra,p,ph
823     :               ,th_wake,tv_wake,lv_wake,cpn_wake
824     :               ,ep,sigp,clw
825     :               ,m,ment,elij,delt,plcl,coef_clos
826     o          ,mp,qp,up,vp,trap,wt,water,evap,b,sigd)
827      endif
828     
829      if (iflag_con.eq.4) then
830       CALL cv_unsat(nloc,ncum,nd,inb,t,q,qs,gz,u,v,p,ph
831     :                   ,h,lv,ep,sigp,clw,m,ment,elij
832     o                   ,iflag,mp,qp,up,vp,wt,water,evap)
833      endif
834c
835      if (debut) THEN
836       print *,'cv_unsat-> '
837       debut=.FALSE.
838      endif !(debut) THEN
839!
840c      print *,'cv_unsat-> mp ',mp
841c      print *,'cv_unsat-> water ',water
842!-------------------------------------------------------------------
843! --- YIELD
844!     (tendencies, precipitation, variables of interface with other
845!      processes, etc)
846!-------------------------------------------------------------------
847
848      if (iflag_con.eq.3) then
849 
850       CALL cv3_yield(nloc,ncum,nd,nd,ntra            ! na->nd
851     :                     ,icb,inb,delt
852     :                     ,t,q,t_wake,q_wake,u,v,tra
853     :                     ,gz,p,ph,h,hp,lv,cpn,th
854     :                     ,ep,clw,m,tp,mp,qp,up,vp,trap
855     :                     ,wt,water,evap,b,sigd
856     :                    ,ment,qent,hent,iflag_mix,uent,vent
857     :                    ,nent,elij,traent,sig
858     :                    ,tv,tvp,wghti
859     :                    ,iflag,precip,Vprecip,ft,fq,fu,fv,ftra
860     :                    ,cbmf,upwd,dnwd,dnwd0,ma,mip
861     :                    ,tls,tps,qcondc,wd
862     :                    ,ftd,fqd)
863!      print *,' cv3_yield -> fqd(1) = ',fqd(1,1)
864      endif
865
866      if (iflag_con.eq.4) then
867       CALL cv_yield(nloc,ncum,nd,nk,icb,inb,delt
868     :              ,t,q,t_wake,q_wake,u,v,tra
869     :              ,gz,p,ph,h,hp,lv,cpn,th
870     :              ,ep,clw,frac,m,mp,qp,up,vp
871     :              ,wt,water,evap
872     :              ,ment,qent,uent,vent,nent,elij
873     :              ,tv,tvp
874     o              ,iflag,wd,qprime,tprime
875     o              ,precip,cbmf,ft,fq,fu,fv,Ma,qcondc)
876      endif
877
878!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
879! --- UNCOMPRESS THE FIELDS
880!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
881
882
883      if (iflag_con.eq.3) then
884       CALL cv3a_uncompress(nloc,len,ncum,nd,ntra,idcum
885     :          ,iflag,icb,inb
886     :          ,precip,sig,w0,ptop2
887     :          ,ft,fq,fu,fv,ftra
888     :          ,Ma,mip,Vprecip,upwd,dnwd,dnwd0
889     ;          ,qcondc,wd,cape,cin
890     :          ,tvp
891     :          ,ftd,fqd
892     :          ,Plim1,Plim2,asupmax,supmax0
893     :          ,asupmaxmin
894     o          ,iflag1,kbas1,ktop1
895     o          ,precip1,sig1,w01,ptop21
896     o          ,ft1,fq1,fu1,fv1,ftra1
897     o          ,Ma1,mip1,Vprecip1,upwd1,dnwd1,dnwd01
898     o          ,qcondc1,wd1,cape1,cin1
899     o          ,tvp1
900     o          ,ftd1,fqd1
901     o          ,Plim11,Plim21,asupmax1,supmax01
902     o          ,asupmaxmin1)
903      endif
904
905      if (iflag_con.eq.4) then
906       CALL cv_uncompress(nloc,len,ncum,nd,idcum
907     :          ,iflag
908     :          ,precip,cbmf
909     :          ,ft,fq,fu,fv
910     :          ,Ma,qcondc
911     o          ,iflag1
912     o          ,precip1,cbmf1
913     o          ,ft1,fq1,fu1,fv1
914     o          ,Ma1,qcondc1 )
915      endif
916
917      ENDIF ! ncum>0
918
9199999  continue
920      return
921      end
922
Note: See TracBrowser for help on using the repository browser.