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

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