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

Last change on this file since 1748 was 1742, checked in by idelkadi, 12 years ago

1- Inclusion des developpements de la these de Romain Pilon sur le
lessivage des aerosols :

a/ par les pluies convectives (modifs cv30_routines et cv3_routines pour

sortir les champs nécessaires au calcul off-line ; modif cvltr)

b/ par les pluies stratiformes (modifs phytrac et introduction

lsc_scav).

2- Choix entre plusieurs schemas pour les pluies stratiformes, commande
par iflag_lscav.

3- Quelques corrections dans la convection "Nouvelle Physique" pour
assurer la conservation des traceurs (cv3p1_mixing et cva_driver) (travail
de Robin Locatelli).

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