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

Last change on this file since 938 was 938, checked in by lmdzadmin, 17 years ago

Enleve prints par defaut
IM

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