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

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

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

  ViewVC Help
Powered by ViewVC 1.1.21