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

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

Suite de la bascule vers une physique avec thermiques, nouvelle convection, poche froide ...
LF

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