source: LMDZ5/branches/LMDZ5_AR5/libf/phylmd/cva_driver.F @ 5447

Last change on this file since 5447 was 1518, checked in by idelkadi, 14 years ago

Modifications des routines de convection :

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