source: LMDZ5/trunk/libf/phylmd/cva_driver.F @ 1694

Last change on this file since 1694 was 1652, checked in by Laurent Fairhead, 12 years ago

Un fichier n'avait pas été commité pour les dernières modifications de A.Cozic


One file not commited in the latest modifications from A. Cozic

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