/[lmdze]/trunk/Sources/phylmd/cv_driver.f
ViewVC logotype

Diff of /trunk/Sources/phylmd/cv_driver.f

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 49 by guez, Wed Aug 24 11:43:14 2011 UTC revision 69 by guez, Mon Feb 18 16:33:12 2013 UTC
# Line 1  Line 1 
1  !  module cv_driver_m
 ! $Header: /home/cvsroot/LMDZ4/libf/phylmd/cv_driver.F,v 1.3 2005/04/15 12:36:17 lmdzadmin Exp $  
 !  
       SUBROUTINE cv_driver(len,nd,ndp1,ntra,iflag_con, &  
                          t1,q1,qs1,u1,v1,tra1, &  
                          p1,ph1,iflag1,ft1,fq1,fu1,fv1,ftra1, &  
                          precip1,VPrecip1, &  
                          cbmf1,sig1,w01, &  
                          icb1,inb1, &  
                          delt,Ma1,upwd1,dnwd1,dnwd01,qcondc1,wd1,cape1, &  
                          da1,phi1,mp1)  
 !  
       use dimens_m  
       use dimphy  
       implicit none  
 !  
 !.............................START PROLOGUE............................  
 !  
 ! PARAMETERS:  
 !      Name            Type         Usage            Description  
 !   ----------      ----------     -------  ----------------------------  
 !  
 !      len           Integer        Input        first (i) dimension  
 !      nd            Integer        Input        vertical (k) dimension  
 !      ndp1          Integer        Input        nd + 1  
 !      ntra          Integer        Input        number of tracors  
 !      iflag_con     Integer        Input        version of convect (3/4)  
 !      t1            Real           Input        temperature  
 !      q1            Real           Input        specific hum  
 !      qs1           Real           Input        sat specific hum  
 !      u1            Real           Input        u-wind  
 !      v1            Real           Input        v-wind  
 !      tra1          Real           Input        tracors  
 !      p1            Real           Input        full level pressure  
 !      ph1           Real           Input        half level pressure  
 !      iflag1        Integer        Output       flag for Emanuel conditions  
 !      ft1           Real           Output       temp tend  
 !      fq1           Real           Output       spec hum tend  
 !      fu1           Real           Output       u-wind tend  
 !      fv1           Real           Output       v-wind tend  
 !      ftra1         Real           Output       tracor tend  
 !      precip1       Real           Output       precipitation  
 !      VPrecip1      Real           Output       vertical profile of precipitations  
 !      cbmf1         Real           Output       cloud base mass flux  
 !      sig1          Real           In/Out       section adiabatic updraft  
 !      w01           Real           In/Out       vertical velocity within adiab updraft  
 !      delt          Real           Input        time step  
 !      Ma1           Real           Output       mass flux adiabatic updraft  
 !      upwd1         Real           Output       total upward mass flux (adiab+mixed)  
 !      dnwd1         Real           Output       saturated downward mass flux (mixed)  
 !      dnwd01        Real           Output       unsaturated downward mass flux  
 !      qcondc1       Real           Output       in-cld mixing ratio of condensed water  
 !      wd1           Real           Output       downdraft velocity scale for sfc fluxes  
 !      cape1         Real           Output       CAPE  
 !  
 ! S. Bony, Mar 2002:  
 !     * Several modules corresponding to different physical processes  
 !     * Several versions of convect may be used:  
 !        - iflag_con=3: version lmd  (previously named convect3)  
 !        - iflag_con=4: version 4.3b (vect. version, previously convect1/2)  
 !   + tard:    - iflag_con=5: version lmd with ice (previously named convectg)  
 ! S. Bony, Oct 2002:  
 !     * Vectorization of convect3 (ie version lmd)  
 !  
 !..............................END PROLOGUE.............................  
 !  
 !  
   
       integer len  
       integer nd  
       integer ndp1  
       integer noff  
       integer, intent(in):: iflag_con  
       integer ntra  
       real t1(len,nd)  
       real q1(len,nd)  
       real qs1(len,nd)  
       real u1(len,nd)  
       real v1(len,nd)  
       real p1(len,nd)  
       real ph1(len,ndp1)  
       integer iflag1(len)  
       real ft1(len,nd)  
       real fq1(len,nd)  
       real fu1(len,nd)  
       real fv1(len,nd)  
       real precip1(len)  
       real cbmf1(len)  
       real VPrecip1(len,nd+1)  
       real Ma1(len,nd)  
       real upwd1(len,nd)  
       real dnwd1(len,nd)  
       real dnwd01(len,nd)  
   
       real qcondc1(len,nd)     ! cld  
       real wd1(len)            ! gust  
       real cape1(len)  
   
       real da1(len,nd),phi1(len,nd,nd),mp1(len,nd)  
       real da(len,nd),phi(len,nd,nd),mp(len,nd)  
       real, intent(in):: tra1(len,nd,ntra)  
       real ftra1(len,nd,ntra)  
   
       real, intent(in):: delt  
   
 !-------------------------------------------------------------------  
 ! --- ARGUMENTS  
 !-------------------------------------------------------------------  
 ! --- On input:  
 !  
 !  t:   Array of absolute temperature (K) of dimension ND, with first  
 !       index corresponding to lowest model level. Note that this array  
 !       will be altered by the subroutine if dry convective adjustment  
 !       occurs and if IPBL is not equal to 0.  
 !  
 !  q:   Array of specific humidity (gm/gm) of dimension ND, with first  
 !       index corresponding to lowest model level. Must be defined  
 !       at same grid levels as T. Note that this array will be altered  
 !       if dry convective adjustment occurs and if IPBL is not equal to 0.  
 !  
 !  qs:  Array of saturation specific humidity of dimension ND, with first  
 !       index corresponding to lowest model level. Must be defined  
 !       at same grid levels as T. Note that this array will be altered  
 !       if dry convective adjustment occurs and if IPBL is not equal to 0.  
 !  
 !  u:   Array of zonal wind velocity (m/s) of dimension ND, witth first  
 !       index corresponding with the lowest model level. Defined at  
 !       same levels as T. Note that this array will be altered if  
 !       dry convective adjustment occurs and if IPBL is not equal to 0.  
 !  
 !  v:   Same as u but for meridional velocity.  
 !  
 !  tra: Array of passive tracer mixing ratio, of dimensions (ND,NTRA),  
 !       where NTRA is the number of different tracers. If no  
 !       convective tracer transport is needed, define a dummy  
 !       input array of dimension (ND,1). Tracers are defined at  
 !       same vertical levels as T. Note that this array will be altered  
 !       if dry convective adjustment occurs and if IPBL is not equal to 0.  
 !  
 !  p:   Array of pressure (mb) of dimension ND, with first  
 !       index corresponding to lowest model level. Must be defined  
 !       at same grid levels as T.  
 !  
 !  ph:  Array of pressure (mb) of dimension ND+1, with first index  
 !       corresponding to lowest level. These pressures are defined at  
 !       levels intermediate between those of P, T, Q and QS. The first  
 !       value of PH should be greater than (i.e. at a lower level than)  
 !       the first value of the array P.  
 !  
 !  nl:  The maximum number of levels to which convection can penetrate, plus 1.  
 !       NL MUST be less than or equal to ND-1.  
 !  
 !  delt: The model time step (sec) between calls to CONVECT  
 !  
 !----------------------------------------------------------------------------  
 ! ---   On Output:  
 !  
 !  iflag: An output integer whose value denotes the following:  
 !       VALUE   INTERPRETATION  
 !       -----   --------------  
 !         0     Moist convection occurs.  
 !         1     Moist convection occurs, but a CFL condition  
 !               on the subsidence warming is violated. This  
 !               does not cause the scheme to terminate.  
 !         2     Moist convection, but no precip because ep(inb) lt 0.0001  
 !         3     No moist convection because new cbmf is 0 and old cbmf is 0.  
 !         4     No moist convection; atmosphere is not  
 !               unstable  
 !         6     No moist convection because ihmin le minorig.  
 !         7     No moist convection because unreasonable  
 !               parcel level temperature or specific humidity.  
 !         8     No moist convection: lifted condensation  
 !               level is above the 200 mb level.  
 !         9     No moist convection: cloud base is higher  
 !               then the level NL-1.  
 !  
 !  ft:   Array of temperature tendency (K/s) of dimension ND, defined at same  
 !        grid levels as T, Q, QS and P.  
 !  
 !  fq:   Array of specific humidity tendencies ((gm/gm)/s) of dimension ND,  
 !        defined at same grid levels as T, Q, QS and P.  
 !  
 !  fu:   Array of forcing of zonal velocity (m/s^2) of dimension ND,  
 !        defined at same grid levels as T.  
 !  
 !  fv:   Same as FU, but for forcing of meridional velocity.  
 !  
 !  ftra: Array of forcing of tracer content, in tracer mixing ratio per  
 !        second, defined at same levels as T. Dimensioned (ND,NTRA).  
 !  
 !  precip: Scalar convective precipitation rate (mm/day).  
 !  
 !  VPrecip: Vertical profile of convective precipitation (kg/m2/s).  
 !  
 !  wd:   A convective downdraft velocity scale. For use in surface  
 !        flux parameterizations. See convect.ps file for details.  
 !  
 !  tprime: A convective downdraft temperature perturbation scale (K).  
 !          For use in surface flux parameterizations. See convect.ps  
 !          file for details.  
 !  
 !  qprime: A convective downdraft specific humidity  
 !          perturbation scale (gm/gm).  
 !          For use in surface flux parameterizations. See convect.ps  
 !          file for details.  
 !  
 !  cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST  
 !        BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT  
 !        ITS NEXT CALL. That is, the value of CBMF must be "remembered"  
 !        by the calling program between calls to CONVECT.  
 !  
 !  det:   Array of detrainment mass flux of dimension ND.  
 !  
 !-------------------------------------------------------------------  
 !  
 !  Local arrays  
 !  
   
       integer i,k,n,il,j  
       integer icbmax  
       integer nk1(klon)  
       integer icb1(klon)  
       integer inb1(klon)  
       integer icbs1(klon)  
   
       real plcl1(klon)  
       real tnk1(klon)  
       real qnk1(klon)  
       real gznk1(klon)  
       real pnk1(klon)  
       real qsnk1(klon)  
       real pbase1(klon)  
       real buoybase1(klon)  
   
       real lv1(klon,klev)  
       real cpn1(klon,klev)  
       real tv1(klon,klev)  
       real gz1(klon,klev)  
       real hm1(klon,klev)  
       real h1(klon,klev)  
       real tp1(klon,klev)  
       real tvp1(klon,klev)  
       real clw1(klon,klev)  
       real sig1(klon,klev)  
       real w01(klon,klev)  
       real th1(klon,klev)  
 !  
       integer ncum  
 !  
 ! (local) compressed fields:  
 !  
       integer nloc  
       parameter (nloc=klon) ! pour l'instant  
   
       integer idcum(nloc)  
       integer iflag(nloc),nk(nloc),icb(nloc)  
       integer nent(nloc,klev)  
       integer icbs(nloc)  
       integer inb(nloc), inbis(nloc)  
   
       real cbmf(nloc),plcl(nloc),tnk(nloc),qnk(nloc),gznk(nloc)  
       real t(nloc,klev),q(nloc,klev),qs(nloc,klev)  
       real u(nloc,klev),v(nloc,klev)  
       real gz(nloc,klev),h(nloc,klev),lv(nloc,klev),cpn(nloc,klev)  
       real p(nloc,klev),ph(nloc,klev+1),tv(nloc,klev),tp(nloc,klev)  
       real clw(nloc,klev)  
       real dph(nloc,klev)  
       real pbase(nloc), buoybase(nloc), th(nloc,klev)  
       real tvp(nloc,klev)  
       real sig(nloc,klev), w0(nloc,klev)  
       real hp(nloc,klev), ep(nloc,klev), sigp(nloc,klev)  
       real frac(nloc), buoy(nloc,klev)  
       real cape(nloc)  
       real m(nloc,klev), ment(nloc,klev,klev), qent(nloc,klev,klev)  
       real uent(nloc,klev,klev), vent(nloc,klev,klev)  
       real ments(nloc,klev,klev), qents(nloc,klev,klev)  
       real sij(nloc,klev,klev), elij(nloc,klev,klev)  
       real qp(nloc,klev), up(nloc,klev), vp(nloc,klev)  
       real wt(nloc,klev), water(nloc,klev), evap(nloc,klev)  
       real b(nloc,klev), ft(nloc,klev), fq(nloc,klev)  
       real fu(nloc,klev), fv(nloc,klev)  
       real upwd(nloc,klev), dnwd(nloc,klev), dnwd0(nloc,klev)  
       real Ma(nloc,klev), mike(nloc,klev), tls(nloc,klev)  
       real tps(nloc,klev), qprime(nloc), tprime(nloc)  
       real precip(nloc)  
       real VPrecip(nloc,klev+1)  
       real tra(nloc,klev,ntra), trap(nloc,klev,ntra)  
       real ftra(nloc,klev,ntra), traent(nloc,klev,klev,ntra)  
       real qcondc(nloc,klev)  ! cld  
       real wd(nloc)           ! gust  
   
 !-------------------------------------------------------------------  
 ! --- SET CONSTANTS AND PARAMETERS  
 !-------------------------------------------------------------------  
   
 ! -- set simulation flags:  
 !   (common cvflag)  
   
        CALL cv_flag  
   
 ! -- set thermodynamical constants:  
 !     (common cvthermo)  
   
        CALL cv_thermo(iflag_con)  
   
 ! -- set convect parameters  
 !  
 !     includes microphysical parameters and parameters that  
 !     control the rate of approach to quasi-equilibrium)  
 !     (common cvparam)  
   
       if (iflag_con.eq.3) then  
        CALL cv3_param(nd,delt)  
       endif  
2    
3        if (iflag_con.eq.4) then    implicit none
4         CALL cv_param(nd)  
5        endif  contains
6    
7  !---------------------------------------------------------------------    SUBROUTINE cv_driver(len, nd, ndp1, ntra, t1, q1, qs1, u1, v1, tra1, p1, &
8  ! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS         ph1, iflag1, ft1, fq1, fu1, fv1, ftra1, precip1, VPrecip1, cbmf1, &
9  !---------------------------------------------------------------------         sig1, w01, icb1, inb1, delt, Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, &
10           cape1, da1, phi1, mp1)
11        do 20 k=1,nd  
12          do 10 i=1,len      ! From LMDZ4/libf/phylmd/cv_driver.F, version 1.3 2005/04/15 12:36:17
13           ft1(i,k)=0.0  
14           fq1(i,k)=0.0      ! Main driver for convection
15           fu1(i,k)=0.0  
16           fv1(i,k)=0.0      use clesphys2, only: iflag_con
17           tvp1(i,k)=0.0      use cv3_param_m, only: cv3_param
18           tp1(i,k)=0.0      USE dimphy, ONLY: klev, klon
19           clw1(i,k)=0.0  
20  !ym      ! PARAMETERS:
21           clw(i,k)=0.0      !      Name            Type         Usage            Description
22           gz1(i,k) = 0.      !   ----------      ----------     -------  ----------------------------
23           VPrecip1(i,k) = 0.  
24           Ma1(i,k)=0.0      !      len           Integer        Input        first (i) dimension
25           upwd1(i,k)=0.0      !      nd            Integer        Input        vertical (k) dimension
26           dnwd1(i,k)=0.0      !      ndp1          Integer        Input        nd + 1
27           dnwd01(i,k)=0.0      !      ntra          Integer        Input        number of tracors
28           qcondc1(i,k)=0.0      !      t1            Real           Input        temperature
29   10     continue      !      q1            Real           Input        specific hum
30   20   continue      !      qs1           Real           Input        sat specific hum
31        !      u1            Real           Input        u-wind
32        do 30 j=1,ntra      !      v1            Real           Input        v-wind
33         do 31 k=1,nd      !      tra1          Real           Input        tracors
34          do 32 i=1,len      !      p1            Real           Input        full level pressure
35           ftra1(i,k,j)=0.0      !      ph1           Real           Input        half level pressure
36   32     continue      !      iflag1        Integer        Output       flag for Emanuel conditions
37   31    continue      !      ft1           Real           Output       temp tend
38   30   continue      !      fq1           Real           Output       spec hum tend
39        !      fu1           Real           Output       u-wind tend
40        do 60 i=1,len      !      fv1           Real           Output       v-wind tend
41          precip1(i)=0.0      !      ftra1         Real           Output       tracor tend
42          iflag1(i)=0      !      precip1       Real           Output       precipitation
43          wd1(i)=0.0      !      VPrecip1      Real           Output       vertical profile of precipitations
44          cape1(i)=0.0      !      cbmf1         Real           Output       cloud base mass flux
45          VPrecip1(i,nd+1)=0.0      !      sig1          Real           In/Out       section adiabatic updraft
46   60   continue      !      w01           Real           In/Out       vertical velocity within adiab updraft
47        !      delt          Real           Input        time step
48        if (iflag_con.eq.3) then      !      Ma1           Real           Output       mass flux adiabatic updraft
49          do il=1,len      !      qcondc1       Real           Output       in-cld mixing ratio of condensed water
50           sig1(il,nd)=sig1(il,nd)+1.      !      wd1           Real           Output       downdraft velocity scale for sfc fluxes
51           sig1(il,nd)=amin1(sig1(il,nd),12.1)      !      cape1         Real           Output       CAPE
52          enddo  
53        endif      ! S. Bony, Mar 2002:
54        !     * Several modules corresponding to different physical processes
55  !--------------------------------------------------------------------      !     * Several versions of convect may be used:
56  ! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY      !        - iflag_con=3: version lmd  (previously named convect3)
57  !--------------------------------------------------------------------      !        - iflag_con=4: version 4.3b (vect. version, previously convect1/2)
58        !   + tard:    - iflag_con=5: version lmd with ice (previously named convectg)
59        if (iflag_con.eq.3) then      ! S. Bony, Oct 2002:
60         CALL cv3_prelim(len,nd,ndp1,t1,q1,p1,ph1             &      !     * Vectorization of convect3 (ie version lmd)
61                       ,lv1,cpn1,tv1,gz1,h1,hm1,th1)! nd->na  
62        endif      integer len
63        integer nd
64        if (iflag_con.eq.4) then      integer ndp1
65         CALL cv_prelim(len,nd,ndp1,t1,q1,p1,ph1 &      integer noff
66                       ,lv1,cpn1,tv1,gz1,h1,hm1)      integer, intent(in):: ntra
67        endif      real, intent(in):: t1(len, nd)
68        real q1(len, nd)
69  !--------------------------------------------------------------------      real qs1(len, nd)
70  ! --- CONVECTIVE FEED      real u1(len, nd)
71  !--------------------------------------------------------------------      real v1(len, nd)
72        real p1(len, nd)
73        if (iflag_con.eq.3) then      real ph1(len, ndp1)
74         CALL cv3_feed(len,nd,t1,q1,qs1,p1,ph1,hm1,gz1            &      integer iflag1(len)
75                 ,nk1,icb1,icbmax,iflag1,tnk1,qnk1,gznk1,plcl1) ! nd->na      real ft1(len, nd)
76        endif      real fq1(len, nd)
77        real fu1(len, nd)
78        if (iflag_con.eq.4) then      real fv1(len, nd)
79         CALL cv_feed(len,nd,t1,q1,qs1,p1,hm1,gz1 &      real precip1(len)
80                 ,nk1,icb1,icbmax,iflag1,tnk1,qnk1,gznk1,plcl1)      real cbmf1(len)
81        endif      real VPrecip1(len, nd+1)
82        real Ma1(len, nd)
83  !--------------------------------------------------------------------      real, intent(out):: upwd1(len, nd) ! total upward mass flux (adiab+mixed)
84  ! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part      real, intent(out):: dnwd1(len, nd) ! saturated downward mass flux (mixed)
85  ! (up through ICB for convect4, up through ICB+1 for convect3)      real, intent(out):: dnwd01(len, nd) ! unsaturated downward mass flux
86  !     Calculates the lifted parcel virtual temperature at nk, the  
87  !     actual temperature, and the adiabatic liquid water content.      real qcondc1(len, nd)     ! cld
88  !--------------------------------------------------------------------      real wd1(len)            ! gust
89        real cape1(len)
90        if (iflag_con.eq.3) then  
91         CALL cv3_undilute1(len,nd,t1,q1,qs1,gz1,plcl1,p1,nk1,icb1   &      real da1(len, nd), phi1(len, nd, nd), mp1(len, nd)
92                                ,tp1,tvp1,clw1,icbs1) ! nd->na      real da(len, nd), phi(len, nd, nd), mp(len, nd)
93        endif      real, intent(in):: tra1(len, nd, ntra)
94        real ftra1(len, nd, ntra)
95        if (iflag_con.eq.4) then  
96         CALL cv_undilute1(len,nd,t1,q1,qs1,gz1,p1,nk1,icb1,icbmax &      real, intent(in):: delt
97                                ,tp1,tvp1,clw1)  
98        endif      !-------------------------------------------------------------------
99        ! --- ARGUMENTS
100  !-------------------------------------------------------------------      !-------------------------------------------------------------------
101  ! --- TRIGGERING      ! --- On input:
102  !-------------------------------------------------------------------  
103        !  t:   Array of absolute temperature (K) of dimension ND, with first
104        if (iflag_con.eq.3) then      !       index corresponding to lowest model level. Note that this array
105         CALL cv3_trigger(len,nd,icb1,plcl1,p1,th1,tv1,tvp1       &      !       will be altered by the subroutine if dry convective adjustment
106                         ,pbase1,buoybase1,iflag1,sig1,w01) ! nd->na      !       occurs and if IPBL is not equal to 0.
107        endif  
108        !  q:   Array of specific humidity (gm/gm) of dimension ND, with first
109        if (iflag_con.eq.4) then      !       index corresponding to lowest model level. Must be defined
110         CALL cv_trigger(len,nd,icb1,cbmf1,tv1,tvp1,iflag1)      !       at same grid levels as T. Note that this array will be altered
111        endif      !       if dry convective adjustment occurs and if IPBL is not equal to 0.
112    
113  !=====================================================================      !  qs:  Array of saturation specific humidity of dimension ND, with first
114  ! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY      !       index corresponding to lowest model level. Must be defined
115  !=====================================================================      !       at same grid levels as T. Note that this array will be altered
116        !       if dry convective adjustment occurs and if IPBL is not equal to 0.
117        ncum=0  
118        do 400 i=1,len      !  u:   Array of zonal wind velocity (m/s) of dimension ND, witth first
119          if(iflag1(i).eq.0)then      !       index corresponding with the lowest model level. Defined at
120             ncum=ncum+1      !       same levels as T. Note that this array will be altered if
121             idcum(ncum)=i      !       dry convective adjustment occurs and if IPBL is not equal to 0.
122          endif  
123   400  continue      !  v:   Same as u but for meridional velocity.
124    
125  !       print*,'klon, ncum = ',len,ncum      !  tra: Array of passive tracer mixing ratio, of dimensions (ND, NTRA),
126        !       where NTRA is the number of different tracers. If no
127        IF (ncum.gt.0) THEN      !       convective tracer transport is needed, define a dummy
128        !       input array of dimension (ND, 1). Tracers are defined at
129  !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^      !       same vertical levels as T. Note that this array will be altered
130  ! --- COMPRESS THE FIELDS      !       if dry convective adjustment occurs and if IPBL is not equal to 0.
131  !        (-> vectorization over convective gridpoints)  
132  !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^      !  p:   Array of pressure (mb) of dimension ND, with first
133        !       index corresponding to lowest model level. Must be defined
134        if (iflag_con.eq.3) then      !       at same grid levels as T.
135         CALL cv3_compress( len,nloc,ncum,nd,ntra &  
136            ,iflag1,nk1,icb1,icbs1 &      !  ph:  Array of pressure (mb) of dimension ND+1, with first index
137            ,plcl1,tnk1,qnk1,gznk1,pbase1,buoybase1 &      !       corresponding to lowest level. These pressures are defined at
138            ,t1,q1,qs1,u1,v1,gz1,th1 &      !       levels intermediate between those of P, T, Q and QS. The first
139            ,tra1 &      !       value of PH should be greater than (i.e. at a lower level than)
140            ,h1,lv1,cpn1,p1,ph1,tv1,tp1,tvp1,clw1  &      !       the first value of the array P.
141            ,sig1,w01 &  
142            ,iflag,nk,icb,icbs &      !  nl:  The maximum number of levels to which convection can penetrate, plus 1.
143            ,plcl,tnk,qnk,gznk,pbase,buoybase &      !       NL MUST be less than or equal to ND-1.
144            ,t,q,qs,u,v,gz,th &  
145            ,tra &      !  delt: The model time step (sec) between calls to CONVECT
146            ,h,lv,cpn,p,ph,tv,tp,tvp,clw  &  
147            ,sig,w0  )      !----------------------------------------------------------------------------
148        endif      ! ---   On Output:
149    
150        if (iflag_con.eq.4) then      !  iflag: An output integer whose value denotes the following:
151         CALL cv_compress( len,nloc,ncum,nd &      !       VALUE   INTERPRETATION
152            ,iflag1,nk1,icb1 &      !       -----   --------------
153            ,cbmf1,plcl1,tnk1,qnk1,gznk1 &      !         0     Moist convection occurs.
154            ,t1,q1,qs1,u1,v1,gz1 &      !         1     Moist convection occurs, but a CFL condition
155            ,h1,lv1,cpn1,p1,ph1,tv1,tp1,tvp1,clw1 &      !               on the subsidence warming is violated. This
156            ,iflag,nk,icb &      !               does not cause the scheme to terminate.
157            ,cbmf,plcl,tnk,qnk,gznk &      !         2     Moist convection, but no precip because ep(inb) lt 0.0001
158            ,t,q,qs,u,v,gz,h,lv,cpn,p,ph,tv,tp,tvp,clw  &      !         3     No moist convection because new cbmf is 0 and old cbmf is 0.
159            ,dph )      !         4     No moist convection; atmosphere is not
160        endif      !               unstable
161        !         6     No moist convection because ihmin le minorig.
162  !-------------------------------------------------------------------      !         7     No moist convection because unreasonable
163  ! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :      !               parcel level temperature or specific humidity.
164  ! ---   FIND THE REST OF THE LIFTED PARCEL TEMPERATURES      !         8     No moist convection: lifted condensation
165  ! ---   &      !               level is above the 200 mb level.
166  ! ---   COMPUTE THE PRECIPITATION EFFICIENCIES AND THE      !         9     No moist convection: cloud base is higher
167  ! ---   FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD      !               then the level NL-1.
168  ! ---   &  
169  ! ---   FIND THE LEVEL OF NEUTRAL BUOYANCY      !  ft:   Array of temperature tendency (K/s) of dimension ND, defined at same
170  !-------------------------------------------------------------------      !        grid levels as T, Q, QS and P.
171    
172        if (iflag_con.eq.3) then      !  fq:   Array of specific humidity tendencies ((gm/gm)/s) of dimension ND,
173         CALL cv3_undilute2(nloc,ncum,nd,icb,icbs,nk         &      !        defined at same grid levels as T, Q, QS and P.
174                                ,tnk,qnk,gznk,t,q,qs,gz &  
175                                ,p,h,tv,lv,pbase,buoybase,plcl &      !  fu:   Array of forcing of zonal velocity (m/s^2) of dimension ND,
176                                ,inb,tp,tvp,clw,hp,ep,sigp,buoy) !na->nd      !        defined at same grid levels as T.
177        endif  
178        !  fv:   Same as FU, but for forcing of meridional velocity.
179        if (iflag_con.eq.4) then  
180         CALL cv_undilute2(nloc,ncum,nd,icb,nk &      !  ftra: Array of forcing of tracer content, in tracer mixing ratio per
181                                ,tnk,qnk,gznk,t,q,qs,gz &      !        second, defined at same levels as T. Dimensioned (ND, NTRA).
182                                ,p,dph,h,tv,lv &  
183                     ,inb,inbis,tp,tvp,clw,hp,ep,sigp,frac)      !  precip: Scalar convective precipitation rate (mm/day).
184        endif  
185        !  VPrecip: Vertical profile of convective precipitation (kg/m2/s).
186  !-------------------------------------------------------------------  
187  ! --- CLOSURE      !  wd:   A convective downdraft velocity scale. For use in surface
188  !-------------------------------------------------------------------      !        flux parameterizations. See convect.ps file for details.
189    
190        if (iflag_con.eq.3) then      !  tprime: A convective downdraft temperature perturbation scale (K).
191         CALL cv3_closure(nloc,ncum,nd,icb,inb               &      !          For use in surface flux parameterizations. See convect.ps
192                               ,pbase,p,ph,tv,buoy &      !          file for details.
193                               ,sig,w0,cape,m) ! na->nd  
194        endif      !  qprime: A convective downdraft specific humidity
195        !          perturbation scale (gm/gm).
196        if (iflag_con.eq.4) then      !          For use in surface flux parameterizations. See convect.ps
197         CALL cv_closure(nloc,ncum,nd,nk,icb &      !          file for details.
198                        ,tv,tvp,p,ph,dph,plcl,cpn &  
199                        ,iflag,cbmf)      !  cbmf: The cloud base mass flux ((kg/m**2)/s). THIS SCALAR VALUE MUST
200        endif      !        BE STORED BY THE CALLING PROGRAM AND RETURNED TO CONVECT AT
201        !        ITS NEXT CALL. That is, the value of CBMF must be "remembered"
202  !-------------------------------------------------------------------      !        by the calling program between calls to CONVECT.
203  ! --- MIXING  
204  !-------------------------------------------------------------------      !  det:   Array of detrainment mass flux of dimension ND.
205    
206        if (iflag_con.eq.3) then      !-------------------------------------------------------------------
207         CALL cv3_mixing(nloc,ncum,nd,nd,ntra,icb,nk,inb     &  
208                             ,ph,t,q,qs,u,v,tra,h,lv,qnk &      !  Local arrays
209                             ,hp,tv,tvp,ep,clw,m,sig &  
210         ,ment,qent,uent,vent, nent,sij,elij,ments,qents,traent)! na->nd      integer i, k, n, il, j
211        endif      integer icbmax
212        integer nk1(klon)
213        if (iflag_con.eq.4) then      integer icb1(klon)
214         CALL cv_mixing(nloc,ncum,nd,icb,nk,inb,inbis &      integer inb1(klon)
215                             ,ph,t,q,qs,u,v,h,lv,qnk &      integer icbs1(klon)
216                             ,hp,tv,tvp,ep,clw,cbmf &  
217                             ,m,ment,qent,uent,vent,nent,sij,elij)      real plcl1(klon)
218        endif      real tnk1(klon)
219        real qnk1(klon)
220  !-------------------------------------------------------------------      real gznk1(klon)
221  ! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS      real pnk1(klon)
222  !-------------------------------------------------------------------      real qsnk1(klon)
223        real pbase1(klon)
224        if (iflag_con.eq.3) then      real buoybase1(klon)
225         CALL cv3_unsat(nloc,ncum,nd,nd,ntra,icb,inb     &  
226                       ,t,q,qs,gz,u,v,tra,p,ph &      real lv1(klon, klev)
227                       ,th,tv,lv,cpn,ep,sigp,clw &      real cpn1(klon, klev)
228                       ,m,ment,elij,delt,plcl &      real tv1(klon, klev)
229                  ,mp,qp,up,vp,trap,wt,water,evap,b)! na->nd      real gz1(klon, klev)
230        endif      real hm1(klon, klev)
231        real h1(klon, klev)
232        if (iflag_con.eq.4) then      real tp1(klon, klev)
233         CALL cv_unsat(nloc,ncum,nd,inb,t,q,qs,gz,u,v,p,ph &      real tvp1(klon, klev)
234                           ,h,lv,ep,sigp,clw,m,ment,elij &      real clw1(klon, klev)
235                           ,iflag,mp,qp,up,vp,wt,water,evap)      real sig1(klon, klev)
236        endif      real w01(klon, klev)
237        real th1(klon, klev)
238  !-------------------------------------------------------------------  
239  ! --- YIELD      integer ncum
240  !     (tendencies, precipitation, variables of interface with other  
241  !      processes, etc)      ! (local) compressed fields:
242  !-------------------------------------------------------------------  
243        integer nloc
244        if (iflag_con.eq.3) then      parameter (nloc=klon) ! pour l'instant
245         CALL cv3_yield(nloc,ncum,nd,nd,ntra             &  
246                             ,icb,inb,delt &      integer idcum(nloc)
247                             ,t,q,u,v,tra,gz,p,ph,h,hp,lv,cpn,th &      integer iflag(nloc), nk(nloc), icb(nloc)
248                             ,ep,clw,m,tp,mp,qp,up,vp,trap &      integer nent(nloc, klev)
249                             ,wt,water,evap,b &      integer icbs(nloc)
250                             ,ment,qent,uent,vent,nent,elij,traent,sig &      integer inb(nloc), inbis(nloc)
251                             ,tv,tvp &  
252                             ,iflag,precip,VPrecip,ft,fq,fu,fv,ftra &      real cbmf(nloc), plcl(nloc), tnk(nloc), qnk(nloc), gznk(nloc)
253                             ,upwd,dnwd,dnwd0,ma,mike,tls,tps,qcondc,wd)! na->nd      real t(nloc, klev), q(nloc, klev), qs(nloc, klev)
254        endif      real u(nloc, klev), v(nloc, klev)
255        real gz(nloc, klev), h(nloc, klev), lv(nloc, klev), cpn(nloc, klev)
256        if (iflag_con.eq.4) then      real p(nloc, klev), ph(nloc, klev+1), tv(nloc, klev), tp(nloc, klev)
257         CALL cv_yield(nloc,ncum,nd,nk,icb,inb,delt &      real clw(nloc, klev)
258                      ,t,q,u,v,gz,p,ph,h,hp,lv,cpn &      real dph(nloc, klev)
259                      ,ep,clw,frac,m,mp,qp,up,vp &      real pbase(nloc), buoybase(nloc), th(nloc, klev)
260                      ,wt,water,evap &      real tvp(nloc, klev)
261                      ,ment,qent,uent,vent,nent,elij &      real sig(nloc, klev), w0(nloc, klev)
262                      ,tv,tvp &      real hp(nloc, klev), ep(nloc, klev), sigp(nloc, klev)
263                      ,iflag,wd,qprime,tprime &      real frac(nloc), buoy(nloc, klev)
264                      ,precip,cbmf,ft,fq,fu,fv,Ma,qcondc)      real cape(nloc)
265        endif      real m(nloc, klev), ment(nloc, klev, klev), qent(nloc, klev, klev)
266        real uent(nloc, klev, klev), vent(nloc, klev, klev)
267  !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^      real ments(nloc, klev, klev), qents(nloc, klev, klev)
268  ! --- passive tracers      real sij(nloc, klev, klev), elij(nloc, klev, klev)
269  !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^      real qp(nloc, klev), up(nloc, klev), vp(nloc, klev)
270        real wt(nloc, klev), water(nloc, klev), evap(nloc, klev)
271        if (iflag_con.eq.3) then      real b(nloc, klev), ft(nloc, klev), fq(nloc, klev)
272         CALL cv3_tracer(nloc,len,ncum,nd,nd, &      real fu(nloc, klev), fv(nloc, klev)
273                          ment,sij,da,phi)      real upwd(nloc, klev), dnwd(nloc, klev), dnwd0(nloc, klev)
274        endif      real Ma(nloc, klev), mike(nloc, klev), tls(nloc, klev)
275        real tps(nloc, klev), qprime(nloc), tprime(nloc)
276  !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^      real precip(nloc)
277  ! --- UNCOMPRESS THE FIELDS      real VPrecip(nloc, klev+1)
278  !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^      real tra(nloc, klev, ntra), trap(nloc, klev, ntra)
279  ! set iflag1 =42 for non convective points      real ftra(nloc, klev, ntra), traent(nloc, klev, klev, ntra)
280        do  i=1,len      real qcondc(nloc, klev)  ! cld
281          iflag1(i)=42      real wd(nloc)           ! gust
282        end do  
283  !      !-------------------------------------------------------------------
284        if (iflag_con.eq.3) then      ! --- SET CONSTANTS AND PARAMETERS
285         CALL cv3_uncompress(nloc,len,ncum,nd,ntra,idcum &      !-------------------------------------------------------------------
286                  ,iflag &  
287                  ,precip,VPrecip,sig,w0 &      ! -- set simulation flags:
288                  ,ft,fq,fu,fv,ftra &      !   (common cvflag)
289                  ,inb  &  
290                  ,Ma,upwd,dnwd,dnwd0,qcondc,wd,cape &      CALL cv_flag
291                  ,da,phi,mp &  
292                  ,iflag1 &      ! -- set thermodynamical constants:
293                  ,precip1,VPrecip1,sig1,w01 &      !     (common cvthermo)
294                  ,ft1,fq1,fu1,fv1,ftra1 &  
295                  ,inb1 &      CALL cv_thermo
296                  ,Ma1,upwd1,dnwd1,dnwd01,qcondc1,wd1,cape1  &  
297                  ,da1,phi1,mp1)      ! -- set convect parameters
298        endif  
299        !     includes microphysical parameters and parameters that
300        if (iflag_con.eq.4) then      !     control the rate of approach to quasi-equilibrium)
301         CALL cv_uncompress(nloc,len,ncum,nd,idcum &      !     (common cvparam)
302                  ,iflag &  
303                  ,precip,cbmf &      if (iflag_con.eq.3) then
304                  ,ft,fq,fu,fv &         CALL cv3_param(nd, delt)
305                  ,Ma,qcondc             &      endif
306                  ,iflag1 &  
307                  ,precip1,cbmf1 &      if (iflag_con.eq.4) then
308                  ,ft1,fq1,fu1,fv1 &         CALL cv_param(nd)
309                  ,Ma1,qcondc1 )      endif
       endif  
310    
311        ENDIF ! ncum>0      !---------------------------------------------------------------------
312        ! --- INITIALIZE OUTPUT ARRAYS AND PARAMETERS
313        !---------------------------------------------------------------------
314    
315        do k=1, nd
316           do  i=1, len
317              ft1(i, k)=0.0
318              fq1(i, k)=0.0
319              fu1(i, k)=0.0
320              fv1(i, k)=0.0
321              tvp1(i, k)=0.0
322              tp1(i, k)=0.0
323              clw1(i, k)=0.0
324              !ym
325              clw(i, k)=0.0
326              gz1(i, k) = 0.
327              VPrecip1(i, k) = 0.
328              Ma1(i, k)=0.0
329              upwd1(i, k)=0.0
330              dnwd1(i, k)=0.0
331              dnwd01(i, k)=0.0
332              qcondc1(i, k)=0.0
333           end do
334        end do
335    
336        do  j=1, ntra
337           do  k=1, nd
338              do  i=1, len
339                 ftra1(i, k, j)=0.0
340              end do
341           end do
342        end do
343    
344        do  i=1, len
345           precip1(i)=0.0
346           iflag1(i)=0
347           wd1(i)=0.0
348           cape1(i)=0.0
349           VPrecip1(i, nd+1)=0.0
350        end do
351    
352        if (iflag_con.eq.3) then
353           do il=1, len
354              sig1(il, nd)=sig1(il, nd)+1.
355              sig1(il, nd)=amin1(sig1(il, nd), 12.1)
356           enddo
357        endif
358    
359        !--------------------------------------------------------------------
360        ! --- CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY
361        !--------------------------------------------------------------------
362    
363        if (iflag_con.eq.3) then
364           CALL cv3_prelim(len, nd, ndp1, t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, &
365                h1, hm1, th1)! nd->na
366        endif
367    
368        if (iflag_con.eq.4) then
369           CALL cv_prelim(len, nd, ndp1, t1, q1, p1, ph1 &
370                , lv1, cpn1, tv1, gz1, h1, hm1)
371        endif
372    
373        !--------------------------------------------------------------------
374        ! --- CONVECTIVE FEED
375        !--------------------------------------------------------------------
376    
377        if (iflag_con.eq.3) then
378           CALL cv3_feed(len, nd, t1, q1, qs1, p1, ph1, hm1, gz1            &
379                , nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! nd->na
380        endif
381    
382        if (iflag_con.eq.4) then
383           CALL cv_feed(len, nd, t1, q1, qs1, p1, hm1, gz1 &
384                , nk1, icb1, icbmax, iflag1, tnk1, qnk1, gznk1, plcl1)
385        endif
386    
387        !--------------------------------------------------------------------
388        ! --- UNDILUTE (ADIABATIC) UPDRAFT / 1st part
389        ! (up through ICB for convect4, up through ICB+1 for convect3)
390        !     Calculates the lifted parcel virtual temperature at nk, the
391        !     actual temperature, and the adiabatic liquid water content.
392        !--------------------------------------------------------------------
393    
394        if (iflag_con.eq.3) then
395           CALL cv3_undilute1(len, nd, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1   &
396                , tp1, tvp1, clw1, icbs1) ! nd->na
397        endif
398    
399        if (iflag_con.eq.4) then
400           CALL cv_undilute1(len, nd, t1, q1, qs1, gz1, p1, nk1, icb1, icbmax &
401                , tp1, tvp1, clw1)
402        endif
403    
404        !-------------------------------------------------------------------
405        ! --- TRIGGERING
406        !-------------------------------------------------------------------
407    
408        if (iflag_con.eq.3) then
409           CALL cv3_trigger(len, nd, icb1, plcl1, p1, th1, tv1, tvp1       &
410                , pbase1, buoybase1, iflag1, sig1, w01) ! nd->na
411        endif
412    
413        if (iflag_con.eq.4) then
414           CALL cv_trigger(len, nd, icb1, cbmf1, tv1, tvp1, iflag1)
415        endif
416    
417        !=====================================================================
418        ! --- IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY
419        !=====================================================================
420    
421        ncum=0
422        do  i=1, len
423           if(iflag1(i).eq.0)then
424              ncum=ncum+1
425              idcum(ncum)=i
426           endif
427        end do
428    
429        !       print*, 'klon, ncum = ', len, ncum
430    
431        IF (ncum.gt.0) THEN
432    
433           !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
434           ! --- COMPRESS THE FIELDS
435           !        (-> vectorization over convective gridpoints)
436           !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
437    
438           if (iflag_con.eq.3) then
439              CALL cv3_compress( len, nloc, ncum, nd, ntra &
440                   , iflag1, nk1, icb1, icbs1 &
441                   , plcl1, tnk1, qnk1, gznk1, pbase1, buoybase1 &
442                   , t1, q1, qs1, u1, v1, gz1, th1 &
443                   , tra1 &
444                   , h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1  &
445                   , sig1, w01 &
446                   , iflag, nk, icb, icbs &
447                   , plcl, tnk, qnk, gznk, pbase, buoybase &
448                   , t, q, qs, u, v, gz, th &
449                   , tra &
450                   , h, lv, cpn, p, ph, tv, tp, tvp, clw  &
451                   , sig, w0  )
452           endif
453    
454           if (iflag_con.eq.4) then
455              CALL cv_compress( len, nloc, ncum, nd &
456                   , iflag1, nk1, icb1 &
457                   , cbmf1, plcl1, tnk1, qnk1, gznk1 &
458                   , t1, q1, qs1, u1, v1, gz1 &
459                   , h1, lv1, cpn1, p1, ph1, tv1, tp1, tvp1, clw1 &
460                   , iflag, nk, icb &
461                   , cbmf, plcl, tnk, qnk, gznk &
462                   , t, q, qs, u, v, gz, h, lv, cpn, p, ph, tv, tp, tvp, clw  &
463                   , dph )
464           endif
465    
466           !-------------------------------------------------------------------
467           ! --- UNDILUTE (ADIABATIC) UPDRAFT / second part :
468           ! ---   FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
469           ! ---   &
470           ! ---   COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
471           ! ---   FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
472           ! ---   &
473           ! ---   FIND THE LEVEL OF NEUTRAL BUOYANCY
474           !-------------------------------------------------------------------
475    
476           if (iflag_con.eq.3) then
477              CALL cv3_undilute2(nloc, ncum, nd, icb, icbs, nk         &
478                   , tnk, qnk, gznk, t, q, qs, gz &
479                   , p, h, tv, lv, pbase, buoybase, plcl &
480                   , inb, tp, tvp, clw, hp, ep, sigp, buoy) !na->nd
481           endif
482    
483           if (iflag_con.eq.4) then
484              CALL cv_undilute2(nloc, ncum, nd, icb, nk &
485                   , tnk, qnk, gznk, t, q, qs, gz &
486                   , p, dph, h, tv, lv &
487                   , inb, inbis, tp, tvp, clw, hp, ep, sigp, frac)
488           endif
489    
490           !-------------------------------------------------------------------
491           ! --- CLOSURE
492           !-------------------------------------------------------------------
493    
494           if (iflag_con.eq.3) then
495              CALL cv3_closure(nloc, ncum, nd, icb, inb               &
496                   , pbase, p, ph, tv, buoy &
497                   , sig, w0, cape, m) ! na->nd
498           endif
499    
500           if (iflag_con.eq.4) then
501              CALL cv_closure(nloc, ncum, nd, nk, icb &
502                   , tv, tvp, p, ph, dph, plcl, cpn &
503                   , iflag, cbmf)
504           endif
505    
506           !-------------------------------------------------------------------
507           ! --- MIXING
508           !-------------------------------------------------------------------
509    
510           if (iflag_con.eq.3) then
511              CALL cv3_mixing(nloc, ncum, nd, nd, ntra, icb, nk, inb     &
512                   , ph, t, q, qs, u, v, tra, h, lv, qnk &
513                   , hp, tv, tvp, ep, clw, m, sig &
514                   , ment, qent, uent, vent, nent, sij, elij, ments, qents, traent)! na->nd
515           endif
516    
517           if (iflag_con.eq.4) then
518              CALL cv_mixing(nloc, ncum, nd, icb, nk, inb, inbis &
519                   , ph, t, q, qs, u, v, h, lv, qnk &
520                   , hp, tv, tvp, ep, clw, cbmf &
521                   , m, ment, qent, uent, vent, nent, sij, elij)
522           endif
523    
524           !-------------------------------------------------------------------
525           ! --- UNSATURATED (PRECIPITATING) DOWNDRAFTS
526           !-------------------------------------------------------------------
527    
528           if (iflag_con.eq.3) then
529              CALL cv3_unsat(nloc, ncum, nd, nd, ntra, icb, inb     &
530                   , t, q, qs, gz, u, v, tra, p, ph &
531                   , th, tv, lv, cpn, ep, sigp, clw &
532                   , m, ment, elij, delt, plcl &
533                   , mp, qp, up, vp, trap, wt, water, evap, b)! na->nd
534           endif
535    
536           if (iflag_con.eq.4) then
537              CALL cv_unsat(nloc, ncum, nd, inb, t, q, qs, gz, u, v, p, ph &
538                   , h, lv, ep, sigp, clw, m, ment, elij &
539                   , iflag, mp, qp, up, vp, wt, water, evap)
540           endif
541    
542           !-------------------------------------------------------------------
543           ! --- YIELD
544           !     (tendencies, precipitation, variables of interface with other
545           !      processes, etc)
546           !-------------------------------------------------------------------
547    
548           if (iflag_con.eq.3) then
549              CALL cv3_yield(nloc, ncum, nd, nd, ntra             &
550                   , icb, inb, delt &
551                   , t, q, u, v, tra, gz, p, ph, h, hp, lv, cpn, th &
552                   , ep, clw, m, tp, mp, qp, up, vp, trap &
553                   , wt, water, evap, b &
554                   , ment, qent, uent, vent, nent, elij, traent, sig &
555                   , tv, tvp &
556                   , iflag, precip, VPrecip, ft, fq, fu, fv, ftra &
557                   , upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc, wd)! na->nd
558           endif
559    
560           if (iflag_con.eq.4) then
561              CALL cv_yield(nloc, ncum, nd, nk, icb, inb, delt &
562                   , t, q, u, v, gz, p, ph, h, hp, lv, cpn &
563                   , ep, clw, frac, m, mp, qp, up, vp &
564                   , wt, water, evap &
565                   , ment, qent, uent, vent, nent, elij &
566                   , tv, tvp &
567                   , iflag, wd, qprime, tprime &
568                   , precip, cbmf, ft, fq, fu, fv, Ma, qcondc)
569           endif
570    
571           !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
572           ! --- passive tracers
573           !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
574    
575           if (iflag_con.eq.3) then
576              CALL cv3_tracer(nloc, len, ncum, nd, nd, &
577                   ment, sij, da, phi)
578           endif
579    
580           !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
581           ! --- UNCOMPRESS THE FIELDS
582           !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
583           ! set iflag1 =42 for non convective points
584           do  i=1, len
585              iflag1(i)=42
586           end do
587    
588           if (iflag_con.eq.3) then
589              CALL cv3_uncompress(nloc, len, ncum, nd, ntra, idcum &
590                   , iflag &
591                   , precip, VPrecip, sig, w0 &
592                   , ft, fq, fu, fv, ftra &
593                   , inb  &
594                   , Ma, upwd, dnwd, dnwd0, qcondc, wd, cape &
595                   , da, phi, mp &
596                   , iflag1 &
597                   , precip1, VPrecip1, sig1, w01 &
598                   , ft1, fq1, fu1, fv1, ftra1 &
599                   , inb1 &
600                   , Ma1, upwd1, dnwd1, dnwd01, qcondc1, wd1, cape1  &
601                   , da1, phi1, mp1)
602           endif
603    
604           if (iflag_con.eq.4) then
605              CALL cv_uncompress(nloc, len, ncum, nd, idcum &
606                   , iflag &
607                   , precip, cbmf &
608                   , ft, fq, fu, fv &
609                   , Ma, qcondc             &
610                   , iflag1 &
611                   , precip1, cbmf1 &
612                   , ft1, fq1, fu1, fv1 &
613                   , Ma1, qcondc1 )
614           endif
615        ENDIF ! ncum>0
616    
617  9999  continue    end SUBROUTINE cv_driver
618    
619        return  end module cv_driver_m
       end  

Legend:
Removed from v.49  
changed lines
  Added in v.69

  ViewVC Help
Powered by ViewVC 1.1.21