/[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

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

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

  ViewVC Help
Powered by ViewVC 1.1.21