source: LMDZ4/branches/LMDZ4V5.0-dev/libf/phylmd/cva_driver.F @ 4138

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